xref: /petsc/src/ts/interface/ts.c (revision be1b0d75ef4789c54b0483a7a345f76eb7155d97)
1 
2 #include <petsc-private/tsimpl.h>        /*I "petscts.h"  I*/
3 #include <petscdmshell.h>
4 #include <petscdmda.h>
5 #include <petscviewer.h>
6 #include <petscdraw.h>
7 
8 /* Logging support */
9 PetscClassId  TS_CLASSID, DMTS_CLASSID;
10 PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
11 
12 const char *const TSExactFinalTimeOptions[] = {"STEPOVER","INTERPOLATE","MATCHSTEP","TSExactFinalTimeOption","TS_EXACTFINALTIME_",0};
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "TSSetTypeFromOptions"
16 /*
17   TSSetTypeFromOptions - Sets the type of ts from user options.
18 
19   Collective on TS
20 
21   Input Parameter:
22 . ts - The ts
23 
24   Level: intermediate
25 
26 .keywords: TS, set, options, database, type
27 .seealso: TSSetFromOptions(), TSSetType()
28 */
29 static PetscErrorCode TSSetTypeFromOptions(TS ts)
30 {
31   PetscBool      opt;
32   const char     *defaultType;
33   char           typeName[256];
34   PetscErrorCode ierr;
35 
36   PetscFunctionBegin;
37   if (((PetscObject)ts)->type_name) defaultType = ((PetscObject)ts)->type_name;
38   else defaultType = TSEULER;
39 
40   if (!TSRegisterAllCalled) {ierr = TSRegisterAll();CHKERRQ(ierr);}
41   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
42   if (opt) {
43     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
44   } else {
45     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
46   }
47   PetscFunctionReturn(0);
48 }
49 
50 struct _n_TSMonitorDrawCtx {
51   PetscViewer   viewer;
52   PetscDrawAxis axis;
53   Vec           initialsolution;
54   PetscBool     showinitial;
55   PetscInt      howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
56   PetscBool     showtimestepandtime;
57   int           color;
58 };
59 
60 #undef __FUNCT__
61 #define __FUNCT__ "TSSetFromOptions"
62 /*@
63    TSSetFromOptions - Sets various TS parameters from user options.
64 
65    Collective on TS
66 
67    Input Parameter:
68 .  ts - the TS context obtained from TSCreate()
69 
70    Options Database Keys:
71 +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP
72 .  -ts_max_steps maxsteps - maximum number of time-steps to take
73 .  -ts_final_time time - maximum time to compute to
74 .  -ts_dt dt - initial time step
75 .  -ts_monitor - print information at each timestep
76 .  -ts_monitor_lg_timestep - Monitor timestep size graphically
77 .  -ts_monitor_lg_solution - Monitor solution graphically
78 .  -ts_monitor_lg_error - Monitor error graphically
79 .  -ts_monitor_lg_snes_iterations - Monitor number nonlinear iterations for each timestep graphically
80 .  -ts_monitor_lg_ksp_iterations - Monitor number nonlinear iterations for each timestep graphically
81 .  -ts_monitor_sp_eig - Monitor eigenvalues of linearized operator graphically
82 .  -ts_monitor_draw_solution - Monitor solution graphically
83 .  -ts_monitor_draw_solution_phase - Monitor solution graphically with phase diagram
84 .  -ts_monitor_draw_error - Monitor error graphically
85 .  -ts_monitor_solution_binary <filename> - Save each solution to a binary file
86 -  -ts_monitor_solution_vtk <filename.vts> - Save each time step to a binary file, use filename-%%03D.vts
87 
88    Developer Note: We should unify all the -ts_monitor options in the way that -xxx_view has been unified
89 
90    Level: beginner
91 
92 .keywords: TS, timestep, set, options, database
93 
94 .seealso: TSGetType()
95 @*/
96 PetscErrorCode  TSSetFromOptions(TS ts)
97 {
98   PetscBool              opt,flg;
99   PetscErrorCode         ierr;
100   PetscViewer            monviewer;
101   char                   monfilename[PETSC_MAX_PATH_LEN];
102   SNES                   snes;
103   TSAdapt                adapt;
104   PetscReal              time_step;
105   TSExactFinalTimeOption eftopt;
106   char                   dir[16];
107 
108   PetscFunctionBegin;
109   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
110   ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr);
111   /* Handle TS type options */
112   ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
113 
114   /* Handle generic TS options */
115   ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,NULL);CHKERRQ(ierr);
116   ierr = PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,NULL);CHKERRQ(ierr);
117   ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,NULL);CHKERRQ(ierr);
118   ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&time_step,&flg);CHKERRQ(ierr);
119   if (flg) {
120     ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
121   }
122   ierr = PetscOptionsEnum("-ts_exact_final_time","Option for handling of final time step","TSSetExactFinalTime",TSExactFinalTimeOptions,(PetscEnum)ts->exact_final_time,(PetscEnum*)&eftopt,&flg);CHKERRQ(ierr);
123   if (flg) {ierr = TSSetExactFinalTime(ts,eftopt);CHKERRQ(ierr);}
124   ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,NULL);CHKERRQ(ierr);
125   ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,NULL);CHKERRQ(ierr);
126   ierr = PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,NULL);CHKERRQ(ierr);
127   ierr = PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,NULL);CHKERRQ(ierr);
128   ierr = PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,NULL);CHKERRQ(ierr);
129 
130   /* Monitor options */
131   ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
132   if (flg) {
133     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),monfilename,&monviewer);CHKERRQ(ierr);
134     ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
135   }
136   ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
137   if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);}
138 
139   ierr = PetscOptionsName("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",&opt);CHKERRQ(ierr);
140   if (opt) {
141     TSMonitorLGCtx ctx;
142     PetscInt       howoften = 1;
143 
144     ierr = PetscOptionsInt("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);CHKERRQ(ierr);
145     ierr = TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
146     ierr = TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
147   }
148   ierr = PetscOptionsName("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",&opt);CHKERRQ(ierr);
149   if (opt) {
150     TSMonitorLGCtx ctx;
151     PetscInt       howoften = 1;
152 
153     ierr = PetscOptionsInt("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",howoften,&howoften,NULL);CHKERRQ(ierr);
154     ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
155     ierr = TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
156   }
157   ierr = PetscOptionsName("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",&opt);CHKERRQ(ierr);
158   if (opt) {
159     TSMonitorLGCtx ctx;
160     PetscInt       howoften = 1;
161 
162     ierr = PetscOptionsInt("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",howoften,&howoften,NULL);CHKERRQ(ierr);
163     ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
164     ierr = TSMonitorSet(ts,TSMonitorLGError,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
165   }
166   ierr = PetscOptionsName("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",&opt);CHKERRQ(ierr);
167   if (opt) {
168     TSMonitorLGCtx ctx;
169     PetscInt       howoften = 1;
170 
171     ierr = PetscOptionsInt("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",howoften,&howoften,NULL);CHKERRQ(ierr);
172     ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
173     ierr = TSMonitorSet(ts,TSMonitorLGSNESIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
174   }
175   ierr = PetscOptionsName("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",&opt);CHKERRQ(ierr);
176   if (opt) {
177     TSMonitorLGCtx ctx;
178     PetscInt       howoften = 1;
179 
180     ierr = PetscOptionsInt("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",howoften,&howoften,NULL);CHKERRQ(ierr);
181     ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
182     ierr = TSMonitorSet(ts,TSMonitorLGKSPIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
183   }
184   ierr = PetscOptionsName("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",&opt);CHKERRQ(ierr);
185   if (opt) {
186     TSMonitorSPEigCtx ctx;
187     PetscInt          howoften = 1;
188 
189     ierr = PetscOptionsInt("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",howoften,&howoften,NULL);CHKERRQ(ierr);
190     ierr = TSMonitorSPEigCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
191     ierr = TSMonitorSet(ts,TSMonitorSPEig,ctx,(PetscErrorCode (*)(void**))TSMonitorSPEigCtxDestroy);CHKERRQ(ierr);
192   }
193   opt  = PETSC_FALSE;
194   ierr = PetscOptionsName("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",&opt);CHKERRQ(ierr);
195   if (opt) {
196     TSMonitorDrawCtx ctx;
197     PetscInt         howoften = 1;
198 
199     ierr = PetscOptionsInt("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",howoften,&howoften,NULL);CHKERRQ(ierr);
200     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
201     ierr = TSMonitorSet(ts,TSMonitorDrawSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
202   }
203   opt  = PETSC_FALSE;
204   ierr = PetscOptionsName("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",&opt);CHKERRQ(ierr);
205   if (opt) {
206     TSMonitorDrawCtx ctx;
207     PetscReal        bounds[4];
208     PetscInt         n = 4;
209     PetscDraw        draw;
210 
211     ierr = PetscOptionsRealArray("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",bounds,&n,NULL);CHKERRQ(ierr);
212     if (n != 4) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Must provide bounding box of phase field");
213     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx);CHKERRQ(ierr);
214     ierr = PetscViewerDrawGetDraw(ctx->viewer,0,&draw);CHKERRQ(ierr);
215     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
216     ierr = PetscDrawAxisCreate(draw,&ctx->axis);CHKERRQ(ierr);
217     ierr = PetscDrawAxisSetLimits(ctx->axis,bounds[0],bounds[2],bounds[1],bounds[3]);CHKERRQ(ierr);
218     ierr = PetscDrawAxisSetLabels(ctx->axis,"Phase Diagram","Variable 1","Variable 2");CHKERRQ(ierr);
219     ierr = PetscDrawAxisDraw(ctx->axis);CHKERRQ(ierr);
220     /* ierr = PetscDrawSetCoordinates(draw,bounds[0],bounds[1],bounds[2],bounds[3]);CHKERRQ(ierr); */
221     ierr = TSMonitorSet(ts,TSMonitorDrawSolutionPhase,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
222   }
223   opt  = PETSC_FALSE;
224   ierr = PetscOptionsName("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",&opt);CHKERRQ(ierr);
225   if (opt) {
226     TSMonitorDrawCtx ctx;
227     PetscInt         howoften = 1;
228 
229     ierr = PetscOptionsInt("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",howoften,&howoften,NULL);CHKERRQ(ierr);
230     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
231     ierr = TSMonitorSet(ts,TSMonitorDrawError,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
232   }
233   opt  = PETSC_FALSE;
234   ierr = PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
235   if (flg) {
236     PetscViewer ctx;
237     if (monfilename[0]) {
238       ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)ts),monfilename,FILE_MODE_WRITE,&ctx);CHKERRQ(ierr);
239       ierr = TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
240     } else {
241       ctx = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)ts));
242       ierr = TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))NULL);CHKERRQ(ierr);
243     }
244   }
245   opt  = PETSC_FALSE;
246   ierr = PetscOptionsString("-ts_monitor_solution_vtk","Save each time step to a binary file, use filename-%%03D.vts","TSMonitorSolutionVTK",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
247   if (flg) {
248     const char *ptr,*ptr2;
249     char       *filetemplate;
250     if (!monfilename[0]) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
251     /* Do some cursory validation of the input. */
252     ierr = PetscStrstr(monfilename,"%",(char**)&ptr);CHKERRQ(ierr);
253     if (!ptr) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
254     for (ptr++; ptr && *ptr; ptr++) {
255       ierr = PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);CHKERRQ(ierr);
256       if (!ptr2 && (*ptr < '0' || '9' < *ptr)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03D.vts");
257       if (ptr2) break;
258     }
259     ierr = PetscStrallocpy(monfilename,&filetemplate);CHKERRQ(ierr);
260     ierr = TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);CHKERRQ(ierr);
261   }
262 
263   ierr = PetscOptionsString("-ts_monitor_dmda_ray","Display a ray of the solution","None","y=0",dir,16,&flg);CHKERRQ(ierr);
264   if (flg) {
265     TSMonitorDMDARayCtx *rayctx;
266     int                 ray = 0;
267     DMDADirection       ddir;
268     DM                  da;
269     PetscMPIInt         rank;
270 
271     if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
272     if (dir[0] == 'x') ddir = DMDA_X;
273     else if (dir[0] == 'y') ddir = DMDA_Y;
274     else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
275     sscanf(dir+2,"%d",&ray);
276 
277     ierr = PetscInfo2(((PetscObject)ts),"Displaying DMDA ray %c = %D\n",dir[0],ray);CHKERRQ(ierr);
278     ierr = PetscNew(TSMonitorDMDARayCtx,&rayctx);CHKERRQ(ierr);
279     ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
280     ierr = DMDAGetRay(da,ddir,ray,&rayctx->ray,&rayctx->scatter);CHKERRQ(ierr);
281     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
282     if (!rank) {
283       ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,600,300,&rayctx->viewer);CHKERRQ(ierr);
284     }
285     ierr = TSMonitorSet(ts,TSMonitorDMDARay,rayctx,TSMonitorDMDARayDestroy);CHKERRQ(ierr);
286   }
287 
288   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
289   ierr = TSAdaptSetFromOptions(adapt);CHKERRQ(ierr);
290 
291   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
292   if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
293 
294   /* Handle specific TS options */
295   if (ts->ops->setfromoptions) {
296     ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
297   }
298 
299   /* process any options handlers added with PetscObjectAddOptionsHandler() */
300   ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
301   ierr = PetscOptionsEnd();CHKERRQ(ierr);
302   PetscFunctionReturn(0);
303 }
304 
305 #undef __FUNCT__
306 #undef __FUNCT__
307 #define __FUNCT__ "TSComputeRHSJacobian"
308 /*@
309    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
310       set with TSSetRHSJacobian().
311 
312    Collective on TS and Vec
313 
314    Input Parameters:
315 +  ts - the TS context
316 .  t - current timestep
317 -  U - input vector
318 
319    Output Parameters:
320 +  A - Jacobian matrix
321 .  B - optional preconditioning matrix
322 -  flag - flag indicating matrix structure
323 
324    Notes:
325    Most users should not need to explicitly call this routine, as it
326    is used internally within the nonlinear solvers.
327 
328    See KSPSetOperators() for important information about setting the
329    flag parameter.
330 
331    Level: developer
332 
333 .keywords: SNES, compute, Jacobian, matrix
334 
335 .seealso:  TSSetRHSJacobian(), KSPSetOperators()
336 @*/
337 PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat *A,Mat *B,MatStructure *flg)
338 {
339   PetscErrorCode ierr;
340   PetscObjectState Ustate;
341   DM             dm;
342   DMTS           tsdm;
343   TSRHSJacobian  rhsjacobianfunc;
344   void           *ctx;
345   TSIJacobian    ijacobianfunc;
346 
347   PetscFunctionBegin;
348   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
349   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
350   PetscCheckSameComm(ts,1,U,3);
351   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
352   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
353   ierr = DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);CHKERRQ(ierr);
354   ierr = DMTSGetIJacobian(dm,&ijacobianfunc,NULL);CHKERRQ(ierr);
355   ierr = PetscObjectStateGet((PetscObject)U,&Ustate);CHKERRQ(ierr);
356   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == U && ts->rhsjacobian.Xstate == Ustate))) {
357     *flg = ts->rhsjacobian.mstructure;
358     PetscFunctionReturn(0);
359   }
360 
361   if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
362 
363   if (ts->rhsjacobian.reuse) {
364     ierr = MatShift(*A,-ts->rhsjacobian.shift);CHKERRQ(ierr);
365     ierr = MatScale(*A,1./ts->rhsjacobian.scale);CHKERRQ(ierr);
366     if (*A != *B) {
367       ierr = MatShift(*B,-ts->rhsjacobian.shift);CHKERRQ(ierr);
368       ierr = MatScale(*B,1./ts->rhsjacobian.scale);CHKERRQ(ierr);
369     }
370     ts->rhsjacobian.shift = 0;
371     ts->rhsjacobian.scale = 1.;
372   }
373 
374   if (rhsjacobianfunc) {
375     ierr = PetscLogEventBegin(TS_JacobianEval,ts,U,*A,*B);CHKERRQ(ierr);
376     *flg = DIFFERENT_NONZERO_PATTERN;
377     PetscStackPush("TS user Jacobian function");
378     ierr = (*rhsjacobianfunc)(ts,t,U,A,B,flg,ctx);CHKERRQ(ierr);
379     PetscStackPop;
380     ierr = PetscLogEventEnd(TS_JacobianEval,ts,U,*A,*B);CHKERRQ(ierr);
381     /* make sure user returned a correct Jacobian and preconditioner */
382     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
383     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
384   } else {
385     ierr = MatZeroEntries(*A);CHKERRQ(ierr);
386     if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);}
387     *flg = SAME_NONZERO_PATTERN;
388   }
389   ts->rhsjacobian.time       = t;
390   ts->rhsjacobian.X          = U;
391   ierr                       = PetscObjectStateGet((PetscObject)U,&ts->rhsjacobian.Xstate);CHKERRQ(ierr);
392   ts->rhsjacobian.mstructure = *flg;
393   PetscFunctionReturn(0);
394 }
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "TSComputeRHSFunction"
398 /*@
399    TSComputeRHSFunction - Evaluates the right-hand-side function.
400 
401    Collective on TS and Vec
402 
403    Input Parameters:
404 +  ts - the TS context
405 .  t - current time
406 -  U - state vector
407 
408    Output Parameter:
409 .  y - right hand side
410 
411    Note:
412    Most users should not need to explicitly call this routine, as it
413    is used internally within the nonlinear solvers.
414 
415    Level: developer
416 
417 .keywords: TS, compute
418 
419 .seealso: TSSetRHSFunction(), TSComputeIFunction()
420 @*/
421 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec U,Vec y)
422 {
423   PetscErrorCode ierr;
424   TSRHSFunction  rhsfunction;
425   TSIFunction    ifunction;
426   void           *ctx;
427   DM             dm;
428 
429   PetscFunctionBegin;
430   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
431   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
432   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
433   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
434   ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctx);CHKERRQ(ierr);
435   ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr);
436 
437   if (!rhsfunction && !ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
438 
439   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,y,0);CHKERRQ(ierr);
440   if (rhsfunction) {
441     PetscStackPush("TS user right-hand-side function");
442     ierr = (*rhsfunction)(ts,t,U,y,ctx);CHKERRQ(ierr);
443     PetscStackPop;
444   } else {
445     ierr = VecZeroEntries(y);CHKERRQ(ierr);
446   }
447 
448   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,y,0);CHKERRQ(ierr);
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "TSComputeSolutionFunction"
454 /*@
455    TSComputeSolutionFunction - Evaluates the solution function.
456 
457    Collective on TS and Vec
458 
459    Input Parameters:
460 +  ts - the TS context
461 -  t - current time
462 
463    Output Parameter:
464 .  U - the solution
465 
466    Note:
467    Most users should not need to explicitly call this routine, as it
468    is used internally within the nonlinear solvers.
469 
470    Level: developer
471 
472 .keywords: TS, compute
473 
474 .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
475 @*/
476 PetscErrorCode TSComputeSolutionFunction(TS ts,PetscReal t,Vec U)
477 {
478   PetscErrorCode     ierr;
479   TSSolutionFunction solutionfunction;
480   void               *ctx;
481   DM                 dm;
482 
483   PetscFunctionBegin;
484   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
485   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
486   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
487   ierr = DMTSGetSolutionFunction(dm,&solutionfunction,&ctx);CHKERRQ(ierr);
488 
489   if (solutionfunction) {
490     PetscStackPush("TS user solution function");
491     ierr = (*solutionfunction)(ts,t,U,ctx);CHKERRQ(ierr);
492     PetscStackPop;
493   }
494   PetscFunctionReturn(0);
495 }
496 #undef __FUNCT__
497 #define __FUNCT__ "TSComputeForcingFunction"
498 /*@
499    TSComputeForcingFunction - Evaluates the forcing function.
500 
501    Collective on TS and Vec
502 
503    Input Parameters:
504 +  ts - the TS context
505 -  t - current time
506 
507    Output Parameter:
508 .  U - the function value
509 
510    Note:
511    Most users should not need to explicitly call this routine, as it
512    is used internally within the nonlinear solvers.
513 
514    Level: developer
515 
516 .keywords: TS, compute
517 
518 .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
519 @*/
520 PetscErrorCode TSComputeForcingFunction(TS ts,PetscReal t,Vec U)
521 {
522   PetscErrorCode     ierr, (*forcing)(TS,PetscReal,Vec,void*);
523   void               *ctx;
524   DM                 dm;
525 
526   PetscFunctionBegin;
527   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
528   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
529   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
530   ierr = DMTSGetForcingFunction(dm,&forcing,&ctx);CHKERRQ(ierr);
531 
532   if (forcing) {
533     PetscStackPush("TS user forcing function");
534     ierr = (*forcing)(ts,t,U,ctx);CHKERRQ(ierr);
535     PetscStackPop;
536   }
537   PetscFunctionReturn(0);
538 }
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "TSGetRHSVec_Private"
542 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
543 {
544   Vec            F;
545   PetscErrorCode ierr;
546 
547   PetscFunctionBegin;
548   *Frhs = NULL;
549   ierr  = TSGetIFunction(ts,&F,NULL,NULL);CHKERRQ(ierr);
550   if (!ts->Frhs) {
551     ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr);
552   }
553   *Frhs = ts->Frhs;
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "TSGetRHSMats_Private"
559 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
560 {
561   Mat            A,B;
562   PetscErrorCode ierr;
563 
564   PetscFunctionBegin;
565   ierr = TSGetIJacobian(ts,&A,&B,NULL,NULL);CHKERRQ(ierr);
566   if (Arhs) {
567     if (!ts->Arhs) {
568       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
569     }
570     *Arhs = ts->Arhs;
571   }
572   if (Brhs) {
573     if (!ts->Brhs) {
574       ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
575     }
576     *Brhs = ts->Brhs;
577   }
578   PetscFunctionReturn(0);
579 }
580 
581 #undef __FUNCT__
582 #define __FUNCT__ "TSComputeIFunction"
583 /*@
584    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,U,Udot)=0
585 
586    Collective on TS and Vec
587 
588    Input Parameters:
589 +  ts - the TS context
590 .  t - current time
591 .  U - state vector
592 .  Udot - time derivative of state vector
593 -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
594 
595    Output Parameter:
596 .  Y - right hand side
597 
598    Note:
599    Most users should not need to explicitly call this routine, as it
600    is used internally within the nonlinear solvers.
601 
602    If the user did did not write their equations in implicit form, this
603    function recasts them in implicit form.
604 
605    Level: developer
606 
607 .keywords: TS, compute
608 
609 .seealso: TSSetIFunction(), TSComputeRHSFunction()
610 @*/
611 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec Y,PetscBool imex)
612 {
613   PetscErrorCode ierr;
614   TSIFunction    ifunction;
615   TSRHSFunction  rhsfunction;
616   void           *ctx;
617   DM             dm;
618 
619   PetscFunctionBegin;
620   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
621   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
622   PetscValidHeaderSpecific(Udot,VEC_CLASSID,4);
623   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
624 
625   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
626   ierr = DMTSGetIFunction(dm,&ifunction,&ctx);CHKERRQ(ierr);
627   ierr = DMTSGetRHSFunction(dm,&rhsfunction,NULL);CHKERRQ(ierr);
628 
629   if (!rhsfunction && !ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
630 
631   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Udot,Y);CHKERRQ(ierr);
632   if (ifunction) {
633     PetscStackPush("TS user implicit function");
634     ierr = (*ifunction)(ts,t,U,Udot,Y,ctx);CHKERRQ(ierr);
635     PetscStackPop;
636   }
637   if (imex) {
638     if (!ifunction) {
639       ierr = VecCopy(Udot,Y);CHKERRQ(ierr);
640     }
641   } else if (rhsfunction) {
642     if (ifunction) {
643       Vec Frhs;
644       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
645       ierr = TSComputeRHSFunction(ts,t,U,Frhs);CHKERRQ(ierr);
646       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
647     } else {
648       ierr = TSComputeRHSFunction(ts,t,U,Y);CHKERRQ(ierr);
649       ierr = VecAYPX(Y,-1,Udot);CHKERRQ(ierr);
650     }
651   }
652   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Udot,Y);CHKERRQ(ierr);
653   PetscFunctionReturn(0);
654 }
655 
656 #undef __FUNCT__
657 #define __FUNCT__ "TSComputeIJacobian"
658 /*@
659    TSComputeIJacobian - Evaluates the Jacobian of the DAE
660 
661    Collective on TS and Vec
662 
663    Input
664       Input Parameters:
665 +  ts - the TS context
666 .  t - current timestep
667 .  U - state vector
668 .  Udot - time derivative of state vector
669 .  shift - shift to apply, see note below
670 -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
671 
672    Output Parameters:
673 +  A - Jacobian matrix
674 .  B - optional preconditioning matrix
675 -  flag - flag indicating matrix structure
676 
677    Notes:
678    If F(t,U,Udot)=0 is the DAE, the required Jacobian is
679 
680    dF/dU + shift*dF/dUdot
681 
682    Most users should not need to explicitly call this routine, as it
683    is used internally within the nonlinear solvers.
684 
685    Level: developer
686 
687 .keywords: TS, compute, Jacobian, matrix
688 
689 .seealso:  TSSetIJacobian()
690 @*/
691 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
692 {
693   PetscErrorCode ierr;
694   TSIJacobian    ijacobian;
695   TSRHSJacobian  rhsjacobian;
696   DM             dm;
697   void           *ctx;
698 
699   PetscFunctionBegin;
700   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
701   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
702   PetscValidHeaderSpecific(Udot,VEC_CLASSID,4);
703   PetscValidPointer(A,6);
704   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
705   PetscValidPointer(B,7);
706   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
707   PetscValidPointer(flg,8);
708 
709   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
710   ierr = DMTSGetIJacobian(dm,&ijacobian,&ctx);CHKERRQ(ierr);
711   ierr = DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);CHKERRQ(ierr);
712 
713   if (!rhsjacobian && !ijacobian) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
714 
715   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
716   ierr = PetscLogEventBegin(TS_JacobianEval,ts,U,*A,*B);CHKERRQ(ierr);
717   if (ijacobian) {
718     *flg = DIFFERENT_NONZERO_PATTERN;
719     PetscStackPush("TS user implicit Jacobian");
720     ierr = (*ijacobian)(ts,t,U,Udot,shift,A,B,flg,ctx);CHKERRQ(ierr);
721     PetscStackPop;
722     /* make sure user returned a correct Jacobian and preconditioner */
723     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
724     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
725   }
726   if (imex) {
727     if (!ijacobian) {  /* system was written as Udot = G(t,U) */
728       ierr = MatZeroEntries(*A);CHKERRQ(ierr);
729       ierr = MatShift(*A,shift);CHKERRQ(ierr);
730       if (*A != *B) {
731         ierr = MatZeroEntries(*B);CHKERRQ(ierr);
732         ierr = MatShift(*B,shift);CHKERRQ(ierr);
733       }
734       *flg = SAME_PRECONDITIONER;
735     }
736   } else {
737     Mat Arhs = NULL,Brhs = NULL;
738     MatStructure flg2;
739     if (rhsjacobian) {
740       ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
741       ierr = TSComputeRHSJacobian(ts,t,U,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
742     }
743     if (Arhs == *A) {           /* No IJacobian, so we only have the RHS matrix */
744       ts->rhsjacobian.scale = -1;
745       ts->rhsjacobian.shift = shift;
746       ierr = MatScale(*A,-1);CHKERRQ(ierr);
747       ierr = MatShift(*A,shift);CHKERRQ(ierr);
748       if (*A != *B) {
749         ierr = MatScale(*B,-1);CHKERRQ(ierr);
750         ierr = MatShift(*B,shift);CHKERRQ(ierr);
751       }
752     } else if (Arhs) {          /* Both IJacobian and RHSJacobian */
753       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
754       if (!ijacobian) {         /* No IJacobian provided, but we have a separate RHS matrix */
755         ierr = MatZeroEntries(*A);CHKERRQ(ierr);
756         ierr = MatShift(*A,shift);CHKERRQ(ierr);
757         if (*A != *B) {
758           ierr = MatZeroEntries(*B);CHKERRQ(ierr);
759           ierr = MatShift(*B,shift);CHKERRQ(ierr);
760         }
761       }
762       ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr);
763       if (*A != *B) {
764         ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr);
765       }
766       *flg = PetscMin(*flg,flg2);
767     }
768   }
769 
770   ierr = PetscLogEventEnd(TS_JacobianEval,ts,U,*A,*B);CHKERRQ(ierr);
771   PetscFunctionReturn(0);
772 }
773 
774 #undef __FUNCT__
775 #define __FUNCT__ "TSSetRHSFunction"
776 /*@C
777     TSSetRHSFunction - Sets the routine for evaluating the function,
778     where U_t = G(t,u).
779 
780     Logically Collective on TS
781 
782     Input Parameters:
783 +   ts - the TS context obtained from TSCreate()
784 .   r - vector to put the computed right hand side (or NULL to have it created)
785 .   f - routine for evaluating the right-hand-side function
786 -   ctx - [optional] user-defined context for private data for the
787           function evaluation routine (may be NULL)
788 
789     Calling sequence of func:
790 $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
791 
792 +   t - current timestep
793 .   u - input vector
794 .   F - function vector
795 -   ctx - [optional] user-defined function context
796 
797     Level: beginner
798 
799 .keywords: TS, timestep, set, right-hand-side, function
800 
801 .seealso: TSSetRHSJacobian(), TSSetIJacobian()
802 @*/
803 PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
804 {
805   PetscErrorCode ierr;
806   SNES           snes;
807   Vec            ralloc = NULL;
808   DM             dm;
809 
810   PetscFunctionBegin;
811   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
812   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
813 
814   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
815   ierr = DMTSSetRHSFunction(dm,f,ctx);CHKERRQ(ierr);
816   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
817   if (!r && !ts->dm && ts->vec_sol) {
818     ierr = VecDuplicate(ts->vec_sol,&ralloc);CHKERRQ(ierr);
819     r    = ralloc;
820   }
821   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
822   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
823   PetscFunctionReturn(0);
824 }
825 
826 #undef __FUNCT__
827 #define __FUNCT__ "TSSetSolutionFunction"
828 /*@C
829     TSSetSolutionFunction - Provide a function that computes the solution of the ODE or DAE
830 
831     Logically Collective on TS
832 
833     Input Parameters:
834 +   ts - the TS context obtained from TSCreate()
835 .   f - routine for evaluating the solution
836 -   ctx - [optional] user-defined context for private data for the
837           function evaluation routine (may be NULL)
838 
839     Calling sequence of func:
840 $     func (TS ts,PetscReal t,Vec u,void *ctx);
841 
842 +   t - current timestep
843 .   u - output vector
844 -   ctx - [optional] user-defined function context
845 
846     Notes:
847     This routine is used for testing accuracy of time integration schemes when you already know the solution.
848     If analytic solutions are not known for your system, consider using the Method of Manufactured Solutions to
849     create closed-form solutions with non-physical forcing terms.
850 
851     For low-dimensional problems solved in serial, such as small discrete systems, TSMonitorLGError() can be used to monitor the error history.
852 
853     Level: beginner
854 
855 .keywords: TS, timestep, set, right-hand-side, function
856 
857 .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetForcingFunction()
858 @*/
859 PetscErrorCode  TSSetSolutionFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
860 {
861   PetscErrorCode ierr;
862   DM             dm;
863 
864   PetscFunctionBegin;
865   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
866   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
867   ierr = DMTSSetSolutionFunction(dm,f,ctx);CHKERRQ(ierr);
868   PetscFunctionReturn(0);
869 }
870 
871 #undef __FUNCT__
872 #define __FUNCT__ "TSSetForcingFunction"
873 /*@C
874     TSSetForcingFunction - Provide a function that computes a forcing term for a ODE or PDE
875 
876     Logically Collective on TS
877 
878     Input Parameters:
879 +   ts - the TS context obtained from TSCreate()
880 .   f - routine for evaluating the forcing function
881 -   ctx - [optional] user-defined context for private data for the
882           function evaluation routine (may be NULL)
883 
884     Calling sequence of func:
885 $     func (TS ts,PetscReal t,Vec u,void *ctx);
886 
887 +   t - current timestep
888 .   u - output vector
889 -   ctx - [optional] user-defined function context
890 
891     Notes:
892     This routine is useful for testing accuracy of time integration schemes when using the Method of Manufactured Solutions to
893     create closed-form solutions with a non-physical forcing term.
894 
895     For low-dimensional problems solved in serial, such as small discrete systems, TSMonitorLGError() can be used to monitor the error history.
896 
897     Level: beginner
898 
899 .keywords: TS, timestep, set, right-hand-side, function
900 
901 .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetSolutionFunction()
902 @*/
903 PetscErrorCode  TSSetForcingFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
904 {
905   PetscErrorCode ierr;
906   DM             dm;
907 
908   PetscFunctionBegin;
909   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
910   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
911   ierr = DMTSSetForcingFunction(dm,f,ctx);CHKERRQ(ierr);
912   PetscFunctionReturn(0);
913 }
914 
915 #undef __FUNCT__
916 #define __FUNCT__ "TSSetRHSJacobian"
917 /*@C
918    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
919    where U_t = G(U,t), as well as the location to store the matrix.
920 
921    Logically Collective on TS
922 
923    Input Parameters:
924 +  ts  - the TS context obtained from TSCreate()
925 .  Amat - (approximate) Jacobian matrix
926 .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
927 .  f   - the Jacobian evaluation routine
928 -  ctx - [optional] user-defined context for private data for the
929          Jacobian evaluation routine (may be NULL)
930 
931    Calling sequence of func:
932 $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
933 
934 +  t - current timestep
935 .  u - input vector
936 .  Amat - (approximate) Jacobian matrix
937 .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
938 .  flag - flag indicating information about the preconditioner matrix
939           structure (same as flag in KSPSetOperators())
940 -  ctx - [optional] user-defined context for matrix evaluation routine
941 
942    Notes:
943    See KSPSetOperators() for important information about setting the flag
944    output parameter in the routine func().  Be sure to read this information!
945 
946    The routine func() takes Mat * as the matrix arguments rather than Mat.
947    This allows the matrix evaluation routine to replace A and/or B with a
948    completely new matrix structure (not just different matrix elements)
949    when appropriate, for instance, if the nonzero structure is changing
950    throughout the global iterations.
951 
952    Level: beginner
953 
954 .keywords: TS, timestep, set, right-hand-side, Jacobian
955 
956 .seealso: SNESComputeJacobianDefaultColor(), TSSetRHSFunction(), TSRHSJacobianSetReuse()
957 
958 @*/
959 PetscErrorCode  TSSetRHSJacobian(TS ts,Mat Amat,Mat Pmat,TSRHSJacobian f,void *ctx)
960 {
961   PetscErrorCode ierr;
962   SNES           snes;
963   DM             dm;
964   TSIJacobian    ijacobian;
965 
966   PetscFunctionBegin;
967   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
968   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
969   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
970   if (Amat) PetscCheckSameComm(ts,1,Amat,2);
971   if (Pmat) PetscCheckSameComm(ts,1,Pmat,3);
972 
973   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
974   ierr = DMTSSetRHSJacobian(dm,f,ctx);CHKERRQ(ierr);
975   if (f == TSComputeRHSJacobianConstant) {
976     /* Handle this case automatically for the user; otherwise user should call themselves. */
977     ierr = TSRHSJacobianSetReuse(ts,PETSC_TRUE);CHKERRQ(ierr);
978   }
979   ierr = DMTSGetIJacobian(dm,&ijacobian,NULL);CHKERRQ(ierr);
980   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
981   if (!ijacobian) {
982     ierr = SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);CHKERRQ(ierr);
983   }
984   if (Amat) {
985     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
986     ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
987 
988     ts->Arhs = Amat;
989   }
990   if (Pmat) {
991     ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);
992     ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
993 
994     ts->Brhs = Pmat;
995   }
996   PetscFunctionReturn(0);
997 }
998 
999 
1000 #undef __FUNCT__
1001 #define __FUNCT__ "TSSetIFunction"
1002 /*@C
1003    TSSetIFunction - Set the function to compute F(t,U,U_t) where F() = 0 is the DAE to be solved.
1004 
1005    Logically Collective on TS
1006 
1007    Input Parameters:
1008 +  ts  - the TS context obtained from TSCreate()
1009 .  r   - vector to hold the residual (or NULL to have it created internally)
1010 .  f   - the function evaluation routine
1011 -  ctx - user-defined context for private data for the function evaluation routine (may be NULL)
1012 
1013    Calling sequence of f:
1014 $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
1015 
1016 +  t   - time at step/stage being solved
1017 .  u   - state vector
1018 .  u_t - time derivative of state vector
1019 .  F   - function vector
1020 -  ctx - [optional] user-defined context for matrix evaluation routine
1021 
1022    Important:
1023    The user MUST call either this routine, TSSetRHSFunction().  This routine must be used when not solving an ODE, for example a DAE.
1024 
1025    Level: beginner
1026 
1027 .keywords: TS, timestep, set, DAE, Jacobian
1028 
1029 .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
1030 @*/
1031 PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
1032 {
1033   PetscErrorCode ierr;
1034   SNES           snes;
1035   Vec            resalloc = NULL;
1036   DM             dm;
1037 
1038   PetscFunctionBegin;
1039   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1040   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
1041 
1042   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1043   ierr = DMTSSetIFunction(dm,f,ctx);CHKERRQ(ierr);
1044 
1045   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1046   if (!res && !ts->dm && ts->vec_sol) {
1047     ierr = VecDuplicate(ts->vec_sol,&resalloc);CHKERRQ(ierr);
1048     res  = resalloc;
1049   }
1050   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
1051   ierr = VecDestroy(&resalloc);CHKERRQ(ierr);
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 #undef __FUNCT__
1056 #define __FUNCT__ "TSGetIFunction"
1057 /*@C
1058    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
1059 
1060    Not Collective
1061 
1062    Input Parameter:
1063 .  ts - the TS context
1064 
1065    Output Parameter:
1066 +  r - vector to hold residual (or NULL)
1067 .  func - the function to compute residual (or NULL)
1068 -  ctx - the function context (or NULL)
1069 
1070    Level: advanced
1071 
1072 .keywords: TS, nonlinear, get, function
1073 
1074 .seealso: TSSetIFunction(), SNESGetFunction()
1075 @*/
1076 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
1077 {
1078   PetscErrorCode ierr;
1079   SNES           snes;
1080   DM             dm;
1081 
1082   PetscFunctionBegin;
1083   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1084   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1085   ierr = SNESGetFunction(snes,r,NULL,NULL);CHKERRQ(ierr);
1086   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1087   ierr = DMTSGetIFunction(dm,func,ctx);CHKERRQ(ierr);
1088   PetscFunctionReturn(0);
1089 }
1090 
1091 #undef __FUNCT__
1092 #define __FUNCT__ "TSGetRHSFunction"
1093 /*@C
1094    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
1095 
1096    Not Collective
1097 
1098    Input Parameter:
1099 .  ts - the TS context
1100 
1101    Output Parameter:
1102 +  r - vector to hold computed right hand side (or NULL)
1103 .  func - the function to compute right hand side (or NULL)
1104 -  ctx - the function context (or NULL)
1105 
1106    Level: advanced
1107 
1108 .keywords: TS, nonlinear, get, function
1109 
1110 .seealso: TSSetRhsfunction(), SNESGetFunction()
1111 @*/
1112 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
1113 {
1114   PetscErrorCode ierr;
1115   SNES           snes;
1116   DM             dm;
1117 
1118   PetscFunctionBegin;
1119   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1120   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1121   ierr = SNESGetFunction(snes,r,NULL,NULL);CHKERRQ(ierr);
1122   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1123   ierr = DMTSGetRHSFunction(dm,func,ctx);CHKERRQ(ierr);
1124   PetscFunctionReturn(0);
1125 }
1126 
1127 #undef __FUNCT__
1128 #define __FUNCT__ "TSSetIJacobian"
1129 /*@C
1130    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
1131         you provided with TSSetIFunction().
1132 
1133    Logically Collective on TS
1134 
1135    Input Parameters:
1136 +  ts  - the TS context obtained from TSCreate()
1137 .  Amat - (approximate) Jacobian matrix
1138 .  Pmat - matrix used to compute preconditioner (usually the same as Amat)
1139 .  f   - the Jacobian evaluation routine
1140 -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)
1141 
1142    Calling sequence of f:
1143 $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *Amat,Mat *Pmat,MatStructure *flag,void *ctx);
1144 
1145 +  t    - time at step/stage being solved
1146 .  U    - state vector
1147 .  U_t  - time derivative of state vector
1148 .  a    - shift
1149 .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1150 .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1151 .  flag - flag indicating information about the preconditioner matrix
1152           structure (same as flag in KSPSetOperators())
1153 -  ctx  - [optional] user-defined context for matrix evaluation routine
1154 
1155    Notes:
1156    The matrices Amat and Pmat are exactly the matrices that are used by SNES for the nonlinear solve.
1157 
1158    The matrix dF/dU + a*dF/dU_t you provide turns out to be
1159    the Jacobian of F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
1160    The time integrator internally approximates U_t by W+a*U where the positive "shift"
1161    a and vector W depend on the integration method, step size, and past states. For example with
1162    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
1163    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
1164 
1165    Level: beginner
1166 
1167 .keywords: TS, timestep, DAE, Jacobian
1168 
1169 .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESComputeJacobianDefaultColor(), SNESComputeJacobianDefault()
1170 
1171 @*/
1172 PetscErrorCode  TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)
1173 {
1174   PetscErrorCode ierr;
1175   SNES           snes;
1176   DM             dm;
1177 
1178   PetscFunctionBegin;
1179   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1180   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1181   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
1182   if (Amat) PetscCheckSameComm(ts,1,Amat,2);
1183   if (Pmat) PetscCheckSameComm(ts,1,Pmat,3);
1184 
1185   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1186   ierr = DMTSSetIJacobian(dm,f,ctx);CHKERRQ(ierr);
1187 
1188   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1189   ierr = SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);CHKERRQ(ierr);
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 #undef __FUNCT__
1194 #define __FUNCT__ "TSRHSJacobianSetReuse"
1195 /*@
1196    TSRHSJacobianSetReuse - restore RHS Jacobian before re-evaluating.  Without this flag, TS will change the sign and
1197    shift the RHS Jacobian for a finite-time-step implicit solve, in which case the user function will need to recompute
1198    the entire Jacobian.  The reuse flag must be set if the evaluation function will assume that the matrix entries have
1199    not been changed by the TS.
1200 
1201    Logically Collective
1202 
1203    Input Arguments:
1204 +  ts - TS context obtained from TSCreate()
1205 -  reuse - PETSC_TRUE if the RHS Jacobian
1206 
1207    Level: intermediate
1208 
1209 .seealso: TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
1210 @*/
1211 PetscErrorCode TSRHSJacobianSetReuse(TS ts,PetscBool reuse)
1212 {
1213   PetscFunctionBegin;
1214   ts->rhsjacobian.reuse = reuse;
1215   PetscFunctionReturn(0);
1216 }
1217 
1218 #undef __FUNCT__
1219 #define __FUNCT__ "TSLoad"
1220 /*@C
1221   TSLoad - Loads a KSP that has been stored in binary  with KSPView().
1222 
1223   Collective on PetscViewer
1224 
1225   Input Parameters:
1226 + newdm - the newly loaded TS, this needs to have been created with TSCreate() or
1227            some related function before a call to TSLoad().
1228 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
1229 
1230    Level: intermediate
1231 
1232   Notes:
1233    The type is determined by the data in the file, any type set into the TS before this call is ignored.
1234 
1235   Notes for advanced users:
1236   Most users should not need to know the details of the binary storage
1237   format, since TSLoad() and TSView() completely hide these details.
1238   But for anyone who's interested, the standard binary matrix storage
1239   format is
1240 .vb
1241      has not yet been determined
1242 .ve
1243 
1244 .seealso: PetscViewerBinaryOpen(), TSView(), MatLoad(), VecLoad()
1245 @*/
1246 PetscErrorCode  TSLoad(TS ts, PetscViewer viewer)
1247 {
1248   PetscErrorCode ierr;
1249   PetscBool      isbinary;
1250   PetscInt       classid;
1251   char           type[256];
1252   DMTS           sdm;
1253   DM             dm;
1254 
1255   PetscFunctionBegin;
1256   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1257   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1258   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1259   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1260 
1261   ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
1262   if (classid != TS_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Not TS next in file");
1263   ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
1264   ierr = TSSetType(ts, type);CHKERRQ(ierr);
1265   if (ts->ops->load) {
1266     ierr = (*ts->ops->load)(ts,viewer);CHKERRQ(ierr);
1267   }
1268   ierr = DMCreate(PetscObjectComm((PetscObject)ts),&dm);CHKERRQ(ierr);
1269   ierr = DMLoad(dm,viewer);CHKERRQ(ierr);
1270   ierr = TSSetDM(ts,dm);CHKERRQ(ierr);
1271   ierr = DMCreateGlobalVector(ts->dm,&ts->vec_sol);CHKERRQ(ierr);
1272   ierr = VecLoad(ts->vec_sol,viewer);CHKERRQ(ierr);
1273   ierr = DMGetDMTS(ts->dm,&sdm);CHKERRQ(ierr);
1274   ierr = DMTSLoad(sdm,viewer);CHKERRQ(ierr);
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 #include <petscdraw.h>
1279 #if defined(PETSC_HAVE_AMS)
1280 #include <petscviewerams.h>
1281 #endif
1282 #undef __FUNCT__
1283 #define __FUNCT__ "TSView"
1284 /*@C
1285     TSView - Prints the TS data structure.
1286 
1287     Collective on TS
1288 
1289     Input Parameters:
1290 +   ts - the TS context obtained from TSCreate()
1291 -   viewer - visualization context
1292 
1293     Options Database Key:
1294 .   -ts_view - calls TSView() at end of TSStep()
1295 
1296     Notes:
1297     The available visualization contexts include
1298 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1299 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1300          output where only the first processor opens
1301          the file.  All other processors send their
1302          data to the first processor to print.
1303 
1304     The user can open an alternative visualization context with
1305     PetscViewerASCIIOpen() - output to a specified file.
1306 
1307     Level: beginner
1308 
1309 .keywords: TS, timestep, view
1310 
1311 .seealso: PetscViewerASCIIOpen()
1312 @*/
1313 PetscErrorCode  TSView(TS ts,PetscViewer viewer)
1314 {
1315   PetscErrorCode ierr;
1316   TSType         type;
1317   PetscBool      iascii,isstring,isundials,isbinary,isdraw;
1318   DMTS           sdm;
1319 #if defined(PETSC_HAVE_AMS)
1320   PetscBool      isams;
1321 #endif
1322 
1323   PetscFunctionBegin;
1324   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1325   if (!viewer) {
1326     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);CHKERRQ(ierr);
1327   }
1328   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1329   PetscCheckSameComm(ts,1,viewer,2);
1330 
1331   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1332   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
1333   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1334   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1335 #if defined(PETSC_HAVE_AMS)
1336   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERAMS,&isams);CHKERRQ(ierr);
1337 #endif
1338   if (iascii) {
1339     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer);CHKERRQ(ierr);
1340     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
1341     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
1342     if (ts->problem_type == TS_NONLINEAR) {
1343       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);CHKERRQ(ierr);
1344       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr);
1345     }
1346     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);CHKERRQ(ierr);
1347     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
1348     ierr = DMGetDMTS(ts->dm,&sdm);CHKERRQ(ierr);
1349     ierr = DMTSView(sdm,viewer);CHKERRQ(ierr);
1350     if (ts->ops->view) {
1351       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1352       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1353       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1354     }
1355   } else if (isstring) {
1356     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
1357     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
1358   } else if (isbinary) {
1359     PetscInt    classid = TS_FILE_CLASSID;
1360     MPI_Comm    comm;
1361     PetscMPIInt rank;
1362     char        type[256];
1363 
1364     ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1365     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1366     if (!rank) {
1367       ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
1368       ierr = PetscStrncpy(type,((PetscObject)ts)->type_name,256);CHKERRQ(ierr);
1369       ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
1370     }
1371     if (ts->ops->view) {
1372       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1373     }
1374     ierr = DMView(ts->dm,viewer);CHKERRQ(ierr);
1375     ierr = VecView(ts->vec_sol,viewer);CHKERRQ(ierr);
1376     ierr = DMGetDMTS(ts->dm,&sdm);CHKERRQ(ierr);
1377     ierr = DMTSView(sdm,viewer);CHKERRQ(ierr);
1378   } else if (isdraw) {
1379     PetscDraw draw;
1380     char      str[36];
1381     PetscReal x,y,bottom,h;
1382 
1383     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1384     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
1385     ierr   = PetscStrcpy(str,"TS: ");CHKERRQ(ierr);
1386     ierr   = PetscStrcat(str,((PetscObject)ts)->type_name);CHKERRQ(ierr);
1387     ierr   = PetscDrawBoxedString(draw,x,y,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
1388     bottom = y - h;
1389     ierr   = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
1390     if (ts->ops->view) {
1391       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1392     }
1393     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
1394 #if defined(PETSC_HAVE_AMS)
1395   } else if (isams) {
1396     if (((PetscObject)ts)->amsmem == -1) {
1397       ierr = PetscObjectViewAMS((PetscObject)ts,viewer);CHKERRQ(ierr);
1398       PetscStackCallAMS(AMS_Memory_take_access,(((PetscObject)ts)->amsmem));
1399       PetscStackCallAMS(AMS_Memory_add_field,(((PetscObject)ts)->amsmem,"time step",&ts->steps,1,AMS_INT,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF));
1400       PetscStackCallAMS(AMS_Memory_add_field,(((PetscObject)ts)->amsmem,"time",&ts->ptime,1,AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF));
1401       PetscStackCallAMS(AMS_Memory_grant_access,(((PetscObject)ts)->amsmem));
1402     }
1403     if (ts->ops->view) {
1404       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1405     }
1406 #endif
1407   }
1408 
1409   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1410   ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
1411   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1412   PetscFunctionReturn(0);
1413 }
1414 
1415 
1416 #undef __FUNCT__
1417 #define __FUNCT__ "TSSetApplicationContext"
1418 /*@
1419    TSSetApplicationContext - Sets an optional user-defined context for
1420    the timesteppers.
1421 
1422    Logically Collective on TS
1423 
1424    Input Parameters:
1425 +  ts - the TS context obtained from TSCreate()
1426 -  usrP - optional user context
1427 
1428    Level: intermediate
1429 
1430 .keywords: TS, timestep, set, application, context
1431 
1432 .seealso: TSGetApplicationContext()
1433 @*/
1434 PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
1435 {
1436   PetscFunctionBegin;
1437   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1438   ts->user = usrP;
1439   PetscFunctionReturn(0);
1440 }
1441 
1442 #undef __FUNCT__
1443 #define __FUNCT__ "TSGetApplicationContext"
1444 /*@
1445     TSGetApplicationContext - Gets the user-defined context for the
1446     timestepper.
1447 
1448     Not Collective
1449 
1450     Input Parameter:
1451 .   ts - the TS context obtained from TSCreate()
1452 
1453     Output Parameter:
1454 .   usrP - user context
1455 
1456     Level: intermediate
1457 
1458 .keywords: TS, timestep, get, application, context
1459 
1460 .seealso: TSSetApplicationContext()
1461 @*/
1462 PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
1463 {
1464   PetscFunctionBegin;
1465   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1466   *(void**)usrP = ts->user;
1467   PetscFunctionReturn(0);
1468 }
1469 
1470 #undef __FUNCT__
1471 #define __FUNCT__ "TSGetTimeStepNumber"
1472 /*@
1473    TSGetTimeStepNumber - Gets the number of time steps completed.
1474 
1475    Not Collective
1476 
1477    Input Parameter:
1478 .  ts - the TS context obtained from TSCreate()
1479 
1480    Output Parameter:
1481 .  iter - number of steps completed so far
1482 
1483    Level: intermediate
1484 
1485 .keywords: TS, timestep, get, iteration, number
1486 .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSSetPostStep()
1487 @*/
1488 PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt *iter)
1489 {
1490   PetscFunctionBegin;
1491   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1492   PetscValidIntPointer(iter,2);
1493   *iter = ts->steps;
1494   PetscFunctionReturn(0);
1495 }
1496 
1497 #undef __FUNCT__
1498 #define __FUNCT__ "TSSetInitialTimeStep"
1499 /*@
1500    TSSetInitialTimeStep - Sets the initial timestep to be used,
1501    as well as the initial time.
1502 
1503    Logically Collective on TS
1504 
1505    Input Parameters:
1506 +  ts - the TS context obtained from TSCreate()
1507 .  initial_time - the initial time
1508 -  time_step - the size of the timestep
1509 
1510    Level: intermediate
1511 
1512 .seealso: TSSetTimeStep(), TSGetTimeStep()
1513 
1514 .keywords: TS, set, initial, timestep
1515 @*/
1516 PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
1517 {
1518   PetscErrorCode ierr;
1519 
1520   PetscFunctionBegin;
1521   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1522   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
1523   ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr);
1524   PetscFunctionReturn(0);
1525 }
1526 
1527 #undef __FUNCT__
1528 #define __FUNCT__ "TSSetTimeStep"
1529 /*@
1530    TSSetTimeStep - Allows one to reset the timestep at any time,
1531    useful for simple pseudo-timestepping codes.
1532 
1533    Logically Collective on TS
1534 
1535    Input Parameters:
1536 +  ts - the TS context obtained from TSCreate()
1537 -  time_step - the size of the timestep
1538 
1539    Level: intermediate
1540 
1541 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1542 
1543 .keywords: TS, set, timestep
1544 @*/
1545 PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
1546 {
1547   PetscFunctionBegin;
1548   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1549   PetscValidLogicalCollectiveReal(ts,time_step,2);
1550   ts->time_step      = time_step;
1551   ts->time_step_orig = time_step;
1552   PetscFunctionReturn(0);
1553 }
1554 
1555 #undef __FUNCT__
1556 #define __FUNCT__ "TSSetExactFinalTime"
1557 /*@
1558    TSSetExactFinalTime - Determines whether to adapt the final time step to
1559      match the exact final time, interpolate solution to the exact final time,
1560      or just return at the final time TS computed.
1561 
1562   Logically Collective on TS
1563 
1564    Input Parameter:
1565 +   ts - the time-step context
1566 -   eftopt - exact final time option
1567 
1568    Level: beginner
1569 
1570 .seealso: TSExactFinalTimeOption
1571 @*/
1572 PetscErrorCode  TSSetExactFinalTime(TS ts,TSExactFinalTimeOption eftopt)
1573 {
1574   PetscFunctionBegin;
1575   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1576   PetscValidLogicalCollectiveEnum(ts,eftopt,2);
1577   ts->exact_final_time = eftopt;
1578   PetscFunctionReturn(0);
1579 }
1580 
1581 #undef __FUNCT__
1582 #define __FUNCT__ "TSGetTimeStep"
1583 /*@
1584    TSGetTimeStep - Gets the current timestep size.
1585 
1586    Not Collective
1587 
1588    Input Parameter:
1589 .  ts - the TS context obtained from TSCreate()
1590 
1591    Output Parameter:
1592 .  dt - the current timestep size
1593 
1594    Level: intermediate
1595 
1596 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1597 
1598 .keywords: TS, get, timestep
1599 @*/
1600 PetscErrorCode  TSGetTimeStep(TS ts,PetscReal *dt)
1601 {
1602   PetscFunctionBegin;
1603   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1604   PetscValidRealPointer(dt,2);
1605   *dt = ts->time_step;
1606   PetscFunctionReturn(0);
1607 }
1608 
1609 #undef __FUNCT__
1610 #define __FUNCT__ "TSGetSolution"
1611 /*@
1612    TSGetSolution - Returns the solution at the present timestep. It
1613    is valid to call this routine inside the function that you are evaluating
1614    in order to move to the new timestep. This vector not changed until
1615    the solution at the next timestep has been calculated.
1616 
1617    Not Collective, but Vec returned is parallel if TS is parallel
1618 
1619    Input Parameter:
1620 .  ts - the TS context obtained from TSCreate()
1621 
1622    Output Parameter:
1623 .  v - the vector containing the solution
1624 
1625    Level: intermediate
1626 
1627 .seealso: TSGetTimeStep()
1628 
1629 .keywords: TS, timestep, get, solution
1630 @*/
1631 PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1632 {
1633   PetscFunctionBegin;
1634   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1635   PetscValidPointer(v,2);
1636   *v = ts->vec_sol;
1637   PetscFunctionReturn(0);
1638 }
1639 
1640 /* ----- Routines to initialize and destroy a timestepper ---- */
1641 #undef __FUNCT__
1642 #define __FUNCT__ "TSSetProblemType"
1643 /*@
1644   TSSetProblemType - Sets the type of problem to be solved.
1645 
1646   Not collective
1647 
1648   Input Parameters:
1649 + ts   - The TS
1650 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1651 .vb
1652          U_t - A U = 0      (linear)
1653          U_t - A(t) U = 0   (linear)
1654          F(t,U,U_t) = 0     (nonlinear)
1655 .ve
1656 
1657    Level: beginner
1658 
1659 .keywords: TS, problem type
1660 .seealso: TSSetUp(), TSProblemType, TS
1661 @*/
1662 PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1663 {
1664   PetscErrorCode ierr;
1665 
1666   PetscFunctionBegin;
1667   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1668   ts->problem_type = type;
1669   if (type == TS_LINEAR) {
1670     SNES snes;
1671     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1672     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
1673   }
1674   PetscFunctionReturn(0);
1675 }
1676 
1677 #undef __FUNCT__
1678 #define __FUNCT__ "TSGetProblemType"
1679 /*@C
1680   TSGetProblemType - Gets the type of problem to be solved.
1681 
1682   Not collective
1683 
1684   Input Parameter:
1685 . ts   - The TS
1686 
1687   Output Parameter:
1688 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1689 .vb
1690          M U_t = A U
1691          M(t) U_t = A(t) U
1692          F(t,U,U_t)
1693 .ve
1694 
1695    Level: beginner
1696 
1697 .keywords: TS, problem type
1698 .seealso: TSSetUp(), TSProblemType, TS
1699 @*/
1700 PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1701 {
1702   PetscFunctionBegin;
1703   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1704   PetscValidIntPointer(type,2);
1705   *type = ts->problem_type;
1706   PetscFunctionReturn(0);
1707 }
1708 
1709 #undef __FUNCT__
1710 #define __FUNCT__ "TSSetUp"
1711 /*@
1712    TSSetUp - Sets up the internal data structures for the later use
1713    of a timestepper.
1714 
1715    Collective on TS
1716 
1717    Input Parameter:
1718 .  ts - the TS context obtained from TSCreate()
1719 
1720    Notes:
1721    For basic use of the TS solvers the user need not explicitly call
1722    TSSetUp(), since these actions will automatically occur during
1723    the call to TSStep().  However, if one wishes to control this
1724    phase separately, TSSetUp() should be called after TSCreate()
1725    and optional routines of the form TSSetXXX(), but before TSStep().
1726 
1727    Level: advanced
1728 
1729 .keywords: TS, timestep, setup
1730 
1731 .seealso: TSCreate(), TSStep(), TSDestroy()
1732 @*/
1733 PetscErrorCode  TSSetUp(TS ts)
1734 {
1735   PetscErrorCode ierr;
1736   DM             dm;
1737   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
1738   PetscErrorCode (*jac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
1739   TSIJacobian    ijac;
1740   TSRHSJacobian  rhsjac;
1741 
1742   PetscFunctionBegin;
1743   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1744   if (ts->setupcalled) PetscFunctionReturn(0);
1745 
1746   if (!((PetscObject)ts)->type_name) {
1747     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1748   }
1749 
1750   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1751 
1752   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
1753 
1754   if (ts->rhsjacobian.reuse) {
1755     Mat Amat,Pmat;
1756     SNES snes;
1757     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1758     ierr = SNESGetJacobian(snes,&Amat,&Pmat,NULL,NULL);CHKERRQ(ierr);
1759     /* Matching matrices implies that an IJacobian is NOT set, because if it had been set, the IJacobian's matrix would
1760      * have displaced the RHS matrix */
1761     if (Amat == ts->Arhs) {
1762       ierr = MatDuplicate(ts->Arhs,MAT_DO_NOT_COPY_VALUES,&Amat);CHKERRQ(ierr);
1763       ierr = SNESSetJacobian(snes,Amat,NULL,NULL,NULL);CHKERRQ(ierr);
1764       ierr = MatDestroy(&Amat);CHKERRQ(ierr);
1765     }
1766     if (Pmat == ts->Brhs) {
1767       ierr = MatDuplicate(ts->Brhs,MAT_DO_NOT_COPY_VALUES,&Pmat);CHKERRQ(ierr);
1768       ierr = SNESSetJacobian(snes,NULL,Pmat,NULL,NULL);CHKERRQ(ierr);
1769       ierr = MatDestroy(&Pmat);CHKERRQ(ierr);
1770     }
1771   }
1772 
1773   if (ts->ops->setup) {
1774     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1775   }
1776 
1777   /* in the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
1778    to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
1779    */
1780   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1781   ierr = DMSNESGetFunction(dm,&func,NULL);CHKERRQ(ierr);
1782   if (!func) {
1783     ierr =DMSNESSetFunction(dm,SNESTSFormFunction,ts);CHKERRQ(ierr);
1784   }
1785   /* if the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
1786      Otherwise, the SNES will use coloring internally to form the Jacobian.
1787    */
1788   ierr = DMSNESGetJacobian(dm,&jac,NULL);CHKERRQ(ierr);
1789   ierr = DMTSGetIJacobian(dm,&ijac,NULL);CHKERRQ(ierr);
1790   ierr = DMTSGetRHSJacobian(dm,&rhsjac,NULL);CHKERRQ(ierr);
1791   if (!jac && (ijac || rhsjac)) {
1792     ierr = DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);CHKERRQ(ierr);
1793   }
1794   ts->setupcalled = PETSC_TRUE;
1795   PetscFunctionReturn(0);
1796 }
1797 
1798 #undef __FUNCT__
1799 #define __FUNCT__ "TSReset"
1800 /*@
1801    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1802 
1803    Collective on TS
1804 
1805    Input Parameter:
1806 .  ts - the TS context obtained from TSCreate()
1807 
1808    Level: beginner
1809 
1810 .keywords: TS, timestep, reset
1811 
1812 .seealso: TSCreate(), TSSetup(), TSDestroy()
1813 @*/
1814 PetscErrorCode  TSReset(TS ts)
1815 {
1816   PetscErrorCode ierr;
1817 
1818   PetscFunctionBegin;
1819   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1820   if (ts->ops->reset) {
1821     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1822   }
1823   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
1824 
1825   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
1826   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1827   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
1828   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1829   ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
1830   ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
1831   ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);
1832 
1833   ts->setupcalled = PETSC_FALSE;
1834   PetscFunctionReturn(0);
1835 }
1836 
1837 #undef __FUNCT__
1838 #define __FUNCT__ "TSDestroy"
1839 /*@
1840    TSDestroy - Destroys the timestepper context that was created
1841    with TSCreate().
1842 
1843    Collective on TS
1844 
1845    Input Parameter:
1846 .  ts - the TS context obtained from TSCreate()
1847 
1848    Level: beginner
1849 
1850 .keywords: TS, timestepper, destroy
1851 
1852 .seealso: TSCreate(), TSSetUp(), TSSolve()
1853 @*/
1854 PetscErrorCode  TSDestroy(TS *ts)
1855 {
1856   PetscErrorCode ierr;
1857 
1858   PetscFunctionBegin;
1859   if (!*ts) PetscFunctionReturn(0);
1860   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
1861   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1862 
1863   ierr = TSReset((*ts));CHKERRQ(ierr);
1864 
1865   /* if memory was published with AMS then destroy it */
1866   ierr = PetscObjectAMSViewOff((PetscObject)*ts);CHKERRQ(ierr);
1867   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
1868 
1869   ierr = TSAdaptDestroy(&(*ts)->adapt);CHKERRQ(ierr);
1870   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
1871   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
1872   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
1873 
1874   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1875   PetscFunctionReturn(0);
1876 }
1877 
1878 #undef __FUNCT__
1879 #define __FUNCT__ "TSGetSNES"
1880 /*@
1881    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1882    a TS (timestepper) context. Valid only for nonlinear problems.
1883 
1884    Not Collective, but SNES is parallel if TS is parallel
1885 
1886    Input Parameter:
1887 .  ts - the TS context obtained from TSCreate()
1888 
1889    Output Parameter:
1890 .  snes - the nonlinear solver context
1891 
1892    Notes:
1893    The user can then directly manipulate the SNES context to set various
1894    options, etc.  Likewise, the user can then extract and manipulate the
1895    KSP, KSP, and PC contexts as well.
1896 
1897    TSGetSNES() does not work for integrators that do not use SNES; in
1898    this case TSGetSNES() returns NULL in snes.
1899 
1900    Level: beginner
1901 
1902 .keywords: timestep, get, SNES
1903 @*/
1904 PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1905 {
1906   PetscErrorCode ierr;
1907 
1908   PetscFunctionBegin;
1909   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1910   PetscValidPointer(snes,2);
1911   if (!ts->snes) {
1912     ierr = SNESCreate(PetscObjectComm((PetscObject)ts),&ts->snes);CHKERRQ(ierr);
1913     ierr = SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);CHKERRQ(ierr);
1914     ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->snes);CHKERRQ(ierr);
1915     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
1916     if (ts->dm) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
1917     if (ts->problem_type == TS_LINEAR) {
1918       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
1919     }
1920   }
1921   *snes = ts->snes;
1922   PetscFunctionReturn(0);
1923 }
1924 
1925 #undef __FUNCT__
1926 #define __FUNCT__ "TSSetSNES"
1927 /*@
1928    TSSetSNES - Set the SNES (nonlinear solver) to be used by the timestepping context
1929 
1930    Collective
1931 
1932    Input Parameter:
1933 +  ts - the TS context obtained from TSCreate()
1934 -  snes - the nonlinear solver context
1935 
1936    Notes:
1937    Most users should have the TS created by calling TSGetSNES()
1938 
1939    Level: developer
1940 
1941 .keywords: timestep, set, SNES
1942 @*/
1943 PetscErrorCode TSSetSNES(TS ts,SNES snes)
1944 {
1945   PetscErrorCode ierr;
1946   PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
1947 
1948   PetscFunctionBegin;
1949   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1950   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
1951   ierr = PetscObjectReference((PetscObject)snes);CHKERRQ(ierr);
1952   ierr = SNESDestroy(&ts->snes);CHKERRQ(ierr);
1953 
1954   ts->snes = snes;
1955 
1956   ierr = SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);CHKERRQ(ierr);
1957   ierr = SNESGetJacobian(ts->snes,NULL,NULL,&func,NULL);CHKERRQ(ierr);
1958   if (func == SNESTSFormJacobian) {
1959     ierr = SNESSetJacobian(ts->snes,NULL,NULL,SNESTSFormJacobian,ts);CHKERRQ(ierr);
1960   }
1961   PetscFunctionReturn(0);
1962 }
1963 
1964 #undef __FUNCT__
1965 #define __FUNCT__ "TSGetKSP"
1966 /*@
1967    TSGetKSP - Returns the KSP (linear solver) associated with
1968    a TS (timestepper) context.
1969 
1970    Not Collective, but KSP is parallel if TS is parallel
1971 
1972    Input Parameter:
1973 .  ts - the TS context obtained from TSCreate()
1974 
1975    Output Parameter:
1976 .  ksp - the nonlinear solver context
1977 
1978    Notes:
1979    The user can then directly manipulate the KSP context to set various
1980    options, etc.  Likewise, the user can then extract and manipulate the
1981    KSP and PC contexts as well.
1982 
1983    TSGetKSP() does not work for integrators that do not use KSP;
1984    in this case TSGetKSP() returns NULL in ksp.
1985 
1986    Level: beginner
1987 
1988 .keywords: timestep, get, KSP
1989 @*/
1990 PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1991 {
1992   PetscErrorCode ierr;
1993   SNES           snes;
1994 
1995   PetscFunctionBegin;
1996   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1997   PetscValidPointer(ksp,2);
1998   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1999   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
2000   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2001   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
2002   PetscFunctionReturn(0);
2003 }
2004 
2005 /* ----------- Routines to set solver parameters ---------- */
2006 
2007 #undef __FUNCT__
2008 #define __FUNCT__ "TSGetDuration"
2009 /*@
2010    TSGetDuration - Gets the maximum number of timesteps to use and
2011    maximum time for iteration.
2012 
2013    Not Collective
2014 
2015    Input Parameters:
2016 +  ts       - the TS context obtained from TSCreate()
2017 .  maxsteps - maximum number of iterations to use, or NULL
2018 -  maxtime  - final time to iterate to, or NULL
2019 
2020    Level: intermediate
2021 
2022 .keywords: TS, timestep, get, maximum, iterations, time
2023 @*/
2024 PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
2025 {
2026   PetscFunctionBegin;
2027   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2028   if (maxsteps) {
2029     PetscValidIntPointer(maxsteps,2);
2030     *maxsteps = ts->max_steps;
2031   }
2032   if (maxtime) {
2033     PetscValidScalarPointer(maxtime,3);
2034     *maxtime = ts->max_time;
2035   }
2036   PetscFunctionReturn(0);
2037 }
2038 
2039 #undef __FUNCT__
2040 #define __FUNCT__ "TSSetDuration"
2041 /*@
2042    TSSetDuration - Sets the maximum number of timesteps to use and
2043    maximum time for iteration.
2044 
2045    Logically Collective on TS
2046 
2047    Input Parameters:
2048 +  ts - the TS context obtained from TSCreate()
2049 .  maxsteps - maximum number of iterations to use
2050 -  maxtime - final time to iterate to
2051 
2052    Options Database Keys:
2053 .  -ts_max_steps <maxsteps> - Sets maxsteps
2054 .  -ts_final_time <maxtime> - Sets maxtime
2055 
2056    Notes:
2057    The default maximum number of iterations is 5000. Default time is 5.0
2058 
2059    Level: intermediate
2060 
2061 .keywords: TS, timestep, set, maximum, iterations
2062 
2063 .seealso: TSSetExactFinalTime()
2064 @*/
2065 PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
2066 {
2067   PetscFunctionBegin;
2068   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2069   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
2070   PetscValidLogicalCollectiveReal(ts,maxtime,2);
2071   if (maxsteps >= 0) ts->max_steps = maxsteps;
2072   if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime;
2073   PetscFunctionReturn(0);
2074 }
2075 
2076 #undef __FUNCT__
2077 #define __FUNCT__ "TSSetSolution"
2078 /*@
2079    TSSetSolution - Sets the initial solution vector
2080    for use by the TS routines.
2081 
2082    Logically Collective on TS and Vec
2083 
2084    Input Parameters:
2085 +  ts - the TS context obtained from TSCreate()
2086 -  u - the solution vector
2087 
2088    Level: beginner
2089 
2090 .keywords: TS, timestep, set, solution, initial conditions
2091 @*/
2092 PetscErrorCode  TSSetSolution(TS ts,Vec u)
2093 {
2094   PetscErrorCode ierr;
2095   DM             dm;
2096 
2097   PetscFunctionBegin;
2098   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2099   PetscValidHeaderSpecific(u,VEC_CLASSID,2);
2100   ierr = PetscObjectReference((PetscObject)u);CHKERRQ(ierr);
2101   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
2102 
2103   ts->vec_sol = u;
2104 
2105   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2106   ierr = DMShellSetGlobalVector(dm,u);CHKERRQ(ierr);
2107   PetscFunctionReturn(0);
2108 }
2109 
2110 #undef __FUNCT__
2111 #define __FUNCT__ "TSSetPreStep"
2112 /*@C
2113   TSSetPreStep - Sets the general-purpose function
2114   called once at the beginning of each time step.
2115 
2116   Logically Collective on TS
2117 
2118   Input Parameters:
2119 + ts   - The TS context obtained from TSCreate()
2120 - func - The function
2121 
2122   Calling sequence of func:
2123 . func (TS ts);
2124 
2125   Level: intermediate
2126 
2127   Note:
2128   If a step is rejected, TSStep() will call this routine again before each attempt.
2129   The last completed time step number can be queried using TSGetTimeStepNumber(), the
2130   size of the step being attempted can be obtained using TSGetTimeStep().
2131 
2132 .keywords: TS, timestep
2133 .seealso: TSSetPreStage(), TSSetPostStage(), TSSetPostStep(), TSStep()
2134 @*/
2135 PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
2136 {
2137   PetscFunctionBegin;
2138   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2139   ts->prestep = func;
2140   PetscFunctionReturn(0);
2141 }
2142 
2143 #undef __FUNCT__
2144 #define __FUNCT__ "TSPreStep"
2145 /*@
2146   TSPreStep - Runs the user-defined pre-step function.
2147 
2148   Collective on TS
2149 
2150   Input Parameters:
2151 . ts   - The TS context obtained from TSCreate()
2152 
2153   Notes:
2154   TSPreStep() is typically used within time stepping implementations,
2155   so most users would not generally call this routine themselves.
2156 
2157   Level: developer
2158 
2159 .keywords: TS, timestep
2160 .seealso: TSSetPreStep(), TSPreStage(), TSPostStage(), TSPostStep()
2161 @*/
2162 PetscErrorCode  TSPreStep(TS ts)
2163 {
2164   PetscErrorCode ierr;
2165 
2166   PetscFunctionBegin;
2167   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2168   if (ts->prestep) {
2169     PetscStackCallStandard((*ts->prestep),(ts));
2170   }
2171   PetscFunctionReturn(0);
2172 }
2173 
2174 #undef __FUNCT__
2175 #define __FUNCT__ "TSSetPreStage"
2176 /*@C
2177   TSSetPreStage - Sets the general-purpose function
2178   called once at the beginning of each stage.
2179 
2180   Logically Collective on TS
2181 
2182   Input Parameters:
2183 + ts   - The TS context obtained from TSCreate()
2184 - func - The function
2185 
2186   Calling sequence of func:
2187 . PetscErrorCode func(TS ts, PetscReal stagetime);
2188 
2189   Level: intermediate
2190 
2191   Note:
2192   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
2193   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
2194   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
2195 
2196 .keywords: TS, timestep
2197 .seealso: TSSetPostStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
2198 @*/
2199 PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
2200 {
2201   PetscFunctionBegin;
2202   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2203   ts->prestage = func;
2204   PetscFunctionReturn(0);
2205 }
2206 
2207 #undef __FUNCT__
2208 #define __FUNCT__ "TSSetPostStage"
2209 /*@C
2210   TSSetPostStage - Sets the general-purpose function
2211   called once at the end of each stage.
2212 
2213   Logically Collective on TS
2214 
2215   Input Parameters:
2216 + ts   - The TS context obtained from TSCreate()
2217 - func - The function
2218 
2219   Calling sequence of func:
2220 . PetscErrorCode func(TS ts, PetscReal stagetime, PetscInt stageindex, Vec* Y);
2221 
2222   Level: intermediate
2223 
2224   Note:
2225   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
2226   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
2227   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
2228 
2229 .keywords: TS, timestep
2230 .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
2231 @*/
2232 PetscErrorCode  TSSetPostStage(TS ts, PetscErrorCode (*func)(TS,PetscReal,PetscInt,Vec*))
2233 {
2234   PetscFunctionBegin;
2235   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2236   ts->poststage = func;
2237   PetscFunctionReturn(0);
2238 }
2239 
2240 #undef __FUNCT__
2241 #define __FUNCT__ "TSPreStage"
2242 /*@
2243   TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()
2244 
2245   Collective on TS
2246 
2247   Input Parameters:
2248 . ts          - The TS context obtained from TSCreate()
2249   stagetime   - The absolute time of the current stage
2250 
2251   Notes:
2252   TSPreStage() is typically used within time stepping implementations,
2253   most users would not generally call this routine themselves.
2254 
2255   Level: developer
2256 
2257 .keywords: TS, timestep
2258 .seealso: TSPostStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
2259 @*/
2260 PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
2261 {
2262   PetscErrorCode ierr;
2263 
2264   PetscFunctionBegin;
2265   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2266   if (ts->prestage) {
2267     PetscStackCallStandard((*ts->prestage),(ts,stagetime));
2268   }
2269   PetscFunctionReturn(0);
2270 }
2271 
2272 #undef __FUNCT__
2273 #define __FUNCT__ "TSPostStage"
2274 /*@
2275   TSPostStage - Runs the user-defined post-stage function set using TSSetPostStage()
2276 
2277   Collective on TS
2278 
2279   Input Parameters:
2280 . ts          - The TS context obtained from TSCreate()
2281   stagetime   - The absolute time of the current stage
2282   stageindex  - Stage number
2283   Y           - Array of vectors (of size = total number
2284                 of stages) with the stage solutions
2285 
2286   Notes:
2287   TSPostStage() is typically used within time stepping implementations,
2288   most users would not generally call this routine themselves.
2289 
2290   Level: developer
2291 
2292 .keywords: TS, timestep
2293 .seealso: TSPreStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
2294 @*/
2295 PetscErrorCode  TSPostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
2296 {
2297   PetscErrorCode ierr;
2298 
2299   PetscFunctionBegin;
2300   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2301   if (ts->prestage) {
2302     PetscStackCallStandard((*ts->poststage),(ts,stagetime,stageindex,Y));
2303   }
2304   PetscFunctionReturn(0);
2305 }
2306 
2307 #undef __FUNCT__
2308 #define __FUNCT__ "TSSetPostStep"
2309 /*@C
2310   TSSetPostStep - Sets the general-purpose function
2311   called once at the end of each time step.
2312 
2313   Logically Collective on TS
2314 
2315   Input Parameters:
2316 + ts   - The TS context obtained from TSCreate()
2317 - func - The function
2318 
2319   Calling sequence of func:
2320 $ func (TS ts);
2321 
2322   Level: intermediate
2323 
2324 .keywords: TS, timestep
2325 .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime()
2326 @*/
2327 PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
2328 {
2329   PetscFunctionBegin;
2330   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2331   ts->poststep = func;
2332   PetscFunctionReturn(0);
2333 }
2334 
2335 #undef __FUNCT__
2336 #define __FUNCT__ "TSPostStep"
2337 /*@
2338   TSPostStep - Runs the user-defined post-step function.
2339 
2340   Collective on TS
2341 
2342   Input Parameters:
2343 . ts   - The TS context obtained from TSCreate()
2344 
2345   Notes:
2346   TSPostStep() is typically used within time stepping implementations,
2347   so most users would not generally call this routine themselves.
2348 
2349   Level: developer
2350 
2351 .keywords: TS, timestep
2352 @*/
2353 PetscErrorCode  TSPostStep(TS ts)
2354 {
2355   PetscErrorCode ierr;
2356 
2357   PetscFunctionBegin;
2358   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2359   if (ts->poststep) {
2360     PetscStackCallStandard((*ts->poststep),(ts));
2361   }
2362   PetscFunctionReturn(0);
2363 }
2364 
2365 /* ------------ Routines to set performance monitoring options ----------- */
2366 
2367 #undef __FUNCT__
2368 #define __FUNCT__ "TSMonitorSet"
2369 /*@C
2370    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
2371    timestep to display the iteration's  progress.
2372 
2373    Logically Collective on TS
2374 
2375    Input Parameters:
2376 +  ts - the TS context obtained from TSCreate()
2377 .  monitor - monitoring routine
2378 .  mctx - [optional] user-defined context for private data for the
2379              monitor routine (use NULL if no context is desired)
2380 -  monitordestroy - [optional] routine that frees monitor context
2381           (may be NULL)
2382 
2383    Calling sequence of monitor:
2384 $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
2385 
2386 +    ts - the TS context
2387 .    steps - iteration number (after the final time step the monitor routine is called with a step of -1, this is at the final time which may have
2388                                been interpolated to)
2389 .    time - current time
2390 .    u - current iterate
2391 -    mctx - [optional] monitoring context
2392 
2393    Notes:
2394    This routine adds an additional monitor to the list of monitors that
2395    already has been loaded.
2396 
2397    Fortran notes: Only a single monitor function can be set for each TS object
2398 
2399    Level: intermediate
2400 
2401 .keywords: TS, timestep, set, monitor
2402 
2403 .seealso: TSMonitorDefault(), TSMonitorCancel()
2404 @*/
2405 PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
2406 {
2407   PetscFunctionBegin;
2408   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2409   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
2410   ts->monitor[ts->numbermonitors]          = monitor;
2411   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
2412   ts->monitorcontext[ts->numbermonitors++] = (void*)mctx;
2413   PetscFunctionReturn(0);
2414 }
2415 
2416 #undef __FUNCT__
2417 #define __FUNCT__ "TSMonitorCancel"
2418 /*@C
2419    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
2420 
2421    Logically Collective on TS
2422 
2423    Input Parameters:
2424 .  ts - the TS context obtained from TSCreate()
2425 
2426    Notes:
2427    There is no way to remove a single, specific monitor.
2428 
2429    Level: intermediate
2430 
2431 .keywords: TS, timestep, set, monitor
2432 
2433 .seealso: TSMonitorDefault(), TSMonitorSet()
2434 @*/
2435 PetscErrorCode  TSMonitorCancel(TS ts)
2436 {
2437   PetscErrorCode ierr;
2438   PetscInt       i;
2439 
2440   PetscFunctionBegin;
2441   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2442   for (i=0; i<ts->numbermonitors; i++) {
2443     if (ts->monitordestroy[i]) {
2444       ierr = (*ts->monitordestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
2445     }
2446   }
2447   ts->numbermonitors = 0;
2448   PetscFunctionReturn(0);
2449 }
2450 
2451 #undef __FUNCT__
2452 #define __FUNCT__ "TSMonitorDefault"
2453 /*@
2454    TSMonitorDefault - Sets the Default monitor
2455 
2456    Level: intermediate
2457 
2458 .keywords: TS, set, monitor
2459 
2460 .seealso: TSMonitorDefault(), TSMonitorSet()
2461 @*/
2462 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
2463 {
2464   PetscErrorCode ierr;
2465   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ts));
2466 
2467   PetscFunctionBegin;
2468   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
2469   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr);
2470   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
2471   PetscFunctionReturn(0);
2472 }
2473 
2474 #undef __FUNCT__
2475 #define __FUNCT__ "TSSetRetainStages"
2476 /*@
2477    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
2478 
2479    Logically Collective on TS
2480 
2481    Input Argument:
2482 .  ts - time stepping context
2483 
2484    Output Argument:
2485 .  flg - PETSC_TRUE or PETSC_FALSE
2486 
2487    Level: intermediate
2488 
2489 .keywords: TS, set
2490 
2491 .seealso: TSInterpolate(), TSSetPostStep()
2492 @*/
2493 PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
2494 {
2495   PetscFunctionBegin;
2496   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2497   ts->retain_stages = flg;
2498   PetscFunctionReturn(0);
2499 }
2500 
2501 #undef __FUNCT__
2502 #define __FUNCT__ "TSInterpolate"
2503 /*@
2504    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
2505 
2506    Collective on TS
2507 
2508    Input Argument:
2509 +  ts - time stepping context
2510 -  t - time to interpolate to
2511 
2512    Output Argument:
2513 .  U - state at given time
2514 
2515    Notes:
2516    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
2517 
2518    Level: intermediate
2519 
2520    Developer Notes:
2521    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
2522 
2523 .keywords: TS, set
2524 
2525 .seealso: TSSetRetainStages(), TSSetPostStep()
2526 @*/
2527 PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec U)
2528 {
2529   PetscErrorCode ierr;
2530 
2531   PetscFunctionBegin;
2532   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2533   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
2534   if (t < ts->ptime - ts->time_step_prev || t > ts->ptime) SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Requested time %G not in last time steps [%G,%G]",t,ts->ptime-ts->time_step_prev,ts->ptime);
2535   if (!ts->ops->interpolate) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
2536   ierr = (*ts->ops->interpolate)(ts,t,U);CHKERRQ(ierr);
2537   PetscFunctionReturn(0);
2538 }
2539 
2540 #undef __FUNCT__
2541 #define __FUNCT__ "TSStep"
2542 /*@
2543    TSStep - Steps one time step
2544 
2545    Collective on TS
2546 
2547    Input Parameter:
2548 .  ts - the TS context obtained from TSCreate()
2549 
2550    Level: intermediate
2551 
2552    Notes:
2553    The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
2554    be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.
2555 
2556    This may over-step the final time provided in TSSetDuration() depending on the time-step used. TSSolve() interpolates to exactly the
2557    time provided in TSSetDuration(). One can use TSInterpolate() to determine an interpolated solution within the final timestep.
2558 
2559 .keywords: TS, timestep, solve
2560 
2561 .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSInterpolate()
2562 @*/
2563 PetscErrorCode  TSStep(TS ts)
2564 {
2565   PetscReal      ptime_prev;
2566   PetscErrorCode ierr;
2567 
2568   PetscFunctionBegin;
2569   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2570   ierr = TSSetUp(ts);CHKERRQ(ierr);
2571 
2572   ts->reason = TS_CONVERGED_ITERATING;
2573   ptime_prev = ts->ptime;
2574 
2575   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
2576   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
2577   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
2578 
2579   ts->time_step_prev = ts->ptime - ptime_prev;
2580 
2581   if (ts->reason < 0) {
2582     if (ts->errorifstepfailed) {
2583       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
2584         SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
2585       } else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
2586     }
2587   } else if (!ts->reason) {
2588     if (ts->steps >= ts->max_steps)     ts->reason = TS_CONVERGED_ITS;
2589     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
2590   }
2591   PetscFunctionReturn(0);
2592 }
2593 
2594 #undef __FUNCT__
2595 #define __FUNCT__ "TSEvaluateStep"
2596 /*@
2597    TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
2598 
2599    Collective on TS
2600 
2601    Input Arguments:
2602 +  ts - time stepping context
2603 .  order - desired order of accuracy
2604 -  done - whether the step was evaluated at this order (pass NULL to generate an error if not available)
2605 
2606    Output Arguments:
2607 .  U - state at the end of the current step
2608 
2609    Level: advanced
2610 
2611    Notes:
2612    This function cannot be called until all stages have been evaluated.
2613    It is normally called by adaptive controllers before a step has been accepted and may also be called by the user after TSStep() has returned.
2614 
2615 .seealso: TSStep(), TSAdapt
2616 @*/
2617 PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec U,PetscBool *done)
2618 {
2619   PetscErrorCode ierr;
2620 
2621   PetscFunctionBegin;
2622   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2623   PetscValidType(ts,1);
2624   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
2625   if (!ts->ops->evaluatestep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
2626   ierr = (*ts->ops->evaluatestep)(ts,order,U,done);CHKERRQ(ierr);
2627   PetscFunctionReturn(0);
2628 }
2629 
2630 #undef __FUNCT__
2631 #define __FUNCT__ "TSSolve"
2632 /*@
2633    TSSolve - Steps the requested number of timesteps.
2634 
2635    Collective on TS
2636 
2637    Input Parameter:
2638 +  ts - the TS context obtained from TSCreate()
2639 -  u - the solution vector  (can be null if TSSetSolution() was used, otherwise must contain the initial conditions)
2640 
2641    Level: beginner
2642 
2643    Notes:
2644    The final time returned by this function may be different from the time of the internally
2645    held state accessible by TSGetSolution() and TSGetTime() because the method may have
2646    stepped over the final time.
2647 
2648 .keywords: TS, timestep, solve
2649 
2650 .seealso: TSCreate(), TSSetSolution(), TSStep()
2651 @*/
2652 PetscErrorCode TSSolve(TS ts,Vec u)
2653 {
2654   PetscBool         flg;
2655   PetscViewer       viewer;
2656   Vec               solution;
2657   PetscErrorCode    ierr;
2658   PetscViewerFormat format;
2659 
2660   PetscFunctionBegin;
2661   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2662   if (u) PetscValidHeaderSpecific(u,VEC_CLASSID,2);
2663   if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
2664     PetscValidHeaderSpecific(u,VEC_CLASSID,2);
2665     if (!ts->vec_sol || u == ts->vec_sol) {
2666       ierr = VecDuplicate(u,&solution);CHKERRQ(ierr);
2667       ierr = TSSetSolution(ts,solution);CHKERRQ(ierr);
2668       ierr = VecDestroy(&solution);CHKERRQ(ierr); /* grant ownership */
2669     }
2670     ierr = VecCopy(u,ts->vec_sol);CHKERRQ(ierr);
2671   } else if (u) {
2672     ierr = TSSetSolution(ts,u);CHKERRQ(ierr);
2673   }
2674   ierr = TSSetUp(ts);CHKERRQ(ierr);
2675   /* reset time step and iteration counters */
2676   ts->steps             = 0;
2677   ts->ksp_its           = 0;
2678   ts->snes_its          = 0;
2679   ts->num_snes_failures = 0;
2680   ts->reject            = 0;
2681   ts->reason            = TS_CONVERGED_ITERATING;
2682 
2683   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject)ts)->prefix,"-ts_view_pre",&viewer,&format,&flg);CHKERRQ(ierr);
2684   if (flg && !PetscPreLoadingOn) {
2685     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
2686     ierr = TSView(ts,viewer);CHKERRQ(ierr);
2687     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
2688     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2689   }
2690 
2691   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
2692     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
2693     ierr = VecCopy(ts->vec_sol,u);CHKERRQ(ierr);
2694     ts->solvetime = ts->ptime;
2695   } else {
2696     /* steps the requested number of timesteps. */
2697     if (ts->steps >= ts->max_steps)     ts->reason = TS_CONVERGED_ITS;
2698     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
2699     while (!ts->reason) {
2700       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2701       ierr = TSStep(ts);CHKERRQ(ierr);
2702       ierr = TSPostStep(ts);CHKERRQ(ierr);
2703     }
2704     if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && ts->ptime > ts->max_time) {
2705       ierr = TSInterpolate(ts,ts->max_time,u);CHKERRQ(ierr);
2706       ts->solvetime = ts->max_time;
2707       solution = u;
2708     } else {
2709       if (u) {ierr = VecCopy(ts->vec_sol,u);CHKERRQ(ierr);}
2710       ts->solvetime = ts->ptime;
2711       solution = ts->vec_sol;
2712     }
2713     ierr = TSMonitor(ts,ts->steps,ts->solvetime,solution);CHKERRQ(ierr);
2714   }
2715   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject)ts)->prefix,"-ts_view",&viewer,&format,&flg);CHKERRQ(ierr);
2716   if (flg && !PetscPreLoadingOn) {
2717     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
2718     ierr = TSView(ts,viewer);CHKERRQ(ierr);
2719     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
2720     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2721   }
2722   PetscFunctionReturn(0);
2723 }
2724 
2725 #undef __FUNCT__
2726 #define __FUNCT__ "TSMonitor"
2727 /*@
2728    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
2729 
2730    Collective on TS
2731 
2732    Input Parameters:
2733 +  ts - time stepping context obtained from TSCreate()
2734 .  step - step number that has just completed
2735 .  ptime - model time of the state
2736 -  u - state at the current model time
2737 
2738    Notes:
2739    TSMonitor() is typically used within the time stepping implementations.
2740    Users might call this function when using the TSStep() interface instead of TSSolve().
2741 
2742    Level: advanced
2743 
2744 .keywords: TS, timestep
2745 @*/
2746 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u)
2747 {
2748   PetscErrorCode ierr;
2749   PetscInt       i,n = ts->numbermonitors;
2750 
2751   PetscFunctionBegin;
2752   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2753   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
2754   for (i=0; i<n; i++) {
2755     ierr = (*ts->monitor[i])(ts,step,ptime,u,ts->monitorcontext[i]);CHKERRQ(ierr);
2756   }
2757   PetscFunctionReturn(0);
2758 }
2759 
2760 /* ------------------------------------------------------------------------*/
2761 #undef __FUNCT__
2762 #define __FUNCT__ "TSMonitorLGCtxCreate"
2763 /*@C
2764    TSMonitorLGCtxCreate - Creates a line graph context for use with
2765    TS to monitor the solution process graphically in various ways
2766 
2767    Collective on TS
2768 
2769    Input Parameters:
2770 +  host - the X display to open, or null for the local machine
2771 .  label - the title to put in the title bar
2772 .  x, y - the screen coordinates of the upper left coordinate of the window
2773 .  m, n - the screen width and height in pixels
2774 -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
2775 
2776    Output Parameter:
2777 .  ctx - the context
2778 
2779    Options Database Key:
2780 +  -ts_monitor_lg_timestep - automatically sets line graph monitor
2781 .  -ts_monitor_lg_solution -
2782 .  -ts_monitor_lg_error -
2783 .  -ts_monitor_lg_ksp_iterations -
2784 .  -ts_monitor_lg_snes_iterations -
2785 -  -lg_indicate_data_points <true,false> - indicate the data points (at each time step) on the plot; default is true
2786 
2787    Notes:
2788    Use TSMonitorLGCtxDestroy() to destroy.
2789 
2790    Level: intermediate
2791 
2792 .keywords: TS, monitor, line graph, residual, seealso
2793 
2794 .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()
2795 
2796 @*/
2797 PetscErrorCode  TSMonitorLGCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtx *ctx)
2798 {
2799   PetscDraw      win;
2800   PetscErrorCode ierr;
2801   PetscBool      flg = PETSC_TRUE;
2802 
2803   PetscFunctionBegin;
2804   ierr = PetscNew(struct _n_TSMonitorLGCtx,ctx);CHKERRQ(ierr);
2805   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr);
2806   ierr = PetscDrawSetFromOptions(win);CHKERRQ(ierr);
2807   ierr = PetscDrawLGCreate(win,1,&(*ctx)->lg);CHKERRQ(ierr);
2808   ierr = PetscOptionsGetBool(NULL,"-lg_indicate_data_points",&flg,NULL);CHKERRQ(ierr);
2809   if (flg) {
2810     ierr = PetscDrawLGIndicateDataPoints((*ctx)->lg);CHKERRQ(ierr);
2811   }
2812   ierr = PetscLogObjectParent((PetscObject)(*ctx)->lg,(PetscObject)win);CHKERRQ(ierr);
2813   (*ctx)->howoften = howoften;
2814   PetscFunctionReturn(0);
2815 }
2816 
2817 #undef __FUNCT__
2818 #define __FUNCT__ "TSMonitorLGTimeStep"
2819 PetscErrorCode TSMonitorLGTimeStep(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
2820 {
2821   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
2822   PetscReal      x   = ptime,y;
2823   PetscErrorCode ierr;
2824 
2825   PetscFunctionBegin;
2826   if (!step) {
2827     PetscDrawAxis axis;
2828     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
2829     ierr = PetscDrawAxisSetLabels(axis,"Timestep as function of time","Time","Time step");CHKERRQ(ierr);
2830     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
2831   }
2832   ierr = TSGetTimeStep(ts,&y);CHKERRQ(ierr);
2833   ierr = PetscDrawLGAddPoint(ctx->lg,&x,&y);CHKERRQ(ierr);
2834   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
2835     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
2836   }
2837   PetscFunctionReturn(0);
2838 }
2839 
2840 #undef __FUNCT__
2841 #define __FUNCT__ "TSMonitorLGCtxDestroy"
2842 /*@C
2843    TSMonitorLGCtxDestroy - Destroys a line graph context that was created
2844    with TSMonitorLGCtxCreate().
2845 
2846    Collective on TSMonitorLGCtx
2847 
2848    Input Parameter:
2849 .  ctx - the monitor context
2850 
2851    Level: intermediate
2852 
2853 .keywords: TS, monitor, line graph, destroy
2854 
2855 .seealso: TSMonitorLGCtxCreate(),  TSMonitorSet(), TSMonitorLGTimeStep();
2856 @*/
2857 PetscErrorCode  TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
2858 {
2859   PetscDraw      draw;
2860   PetscErrorCode ierr;
2861 
2862   PetscFunctionBegin;
2863   ierr = PetscDrawLGGetDraw((*ctx)->lg,&draw);CHKERRQ(ierr);
2864   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
2865   ierr = PetscDrawLGDestroy(&(*ctx)->lg);CHKERRQ(ierr);
2866   ierr = PetscFree(*ctx);CHKERRQ(ierr);
2867   PetscFunctionReturn(0);
2868 }
2869 
2870 #undef __FUNCT__
2871 #define __FUNCT__ "TSGetTime"
2872 /*@
2873    TSGetTime - Gets the time of the most recently completed step.
2874 
2875    Not Collective
2876 
2877    Input Parameter:
2878 .  ts - the TS context obtained from TSCreate()
2879 
2880    Output Parameter:
2881 .  t  - the current time
2882 
2883    Level: beginner
2884 
2885    Note:
2886    When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
2887    TSSetPreStage(), TSSetPostStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.
2888 
2889 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
2890 
2891 .keywords: TS, get, time
2892 @*/
2893 PetscErrorCode  TSGetTime(TS ts,PetscReal *t)
2894 {
2895   PetscFunctionBegin;
2896   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2897   PetscValidRealPointer(t,2);
2898   *t = ts->ptime;
2899   PetscFunctionReturn(0);
2900 }
2901 
2902 #undef __FUNCT__
2903 #define __FUNCT__ "TSSetTime"
2904 /*@
2905    TSSetTime - Allows one to reset the time.
2906 
2907    Logically Collective on TS
2908 
2909    Input Parameters:
2910 +  ts - the TS context obtained from TSCreate()
2911 -  time - the time
2912 
2913    Level: intermediate
2914 
2915 .seealso: TSGetTime(), TSSetDuration()
2916 
2917 .keywords: TS, set, time
2918 @*/
2919 PetscErrorCode  TSSetTime(TS ts, PetscReal t)
2920 {
2921   PetscFunctionBegin;
2922   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2923   PetscValidLogicalCollectiveReal(ts,t,2);
2924   ts->ptime = t;
2925   PetscFunctionReturn(0);
2926 }
2927 
2928 #undef __FUNCT__
2929 #define __FUNCT__ "TSSetOptionsPrefix"
2930 /*@C
2931    TSSetOptionsPrefix - Sets the prefix used for searching for all
2932    TS options in the database.
2933 
2934    Logically Collective on TS
2935 
2936    Input Parameter:
2937 +  ts     - The TS context
2938 -  prefix - The prefix to prepend to all option names
2939 
2940    Notes:
2941    A hyphen (-) must NOT be given at the beginning of the prefix name.
2942    The first character of all runtime options is AUTOMATICALLY the
2943    hyphen.
2944 
2945    Level: advanced
2946 
2947 .keywords: TS, set, options, prefix, database
2948 
2949 .seealso: TSSetFromOptions()
2950 
2951 @*/
2952 PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2953 {
2954   PetscErrorCode ierr;
2955   SNES           snes;
2956 
2957   PetscFunctionBegin;
2958   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2959   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2960   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2961   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2962   PetscFunctionReturn(0);
2963 }
2964 
2965 
2966 #undef __FUNCT__
2967 #define __FUNCT__ "TSAppendOptionsPrefix"
2968 /*@C
2969    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2970    TS options in the database.
2971 
2972    Logically Collective on TS
2973 
2974    Input Parameter:
2975 +  ts     - The TS context
2976 -  prefix - The prefix to prepend to all option names
2977 
2978    Notes:
2979    A hyphen (-) must NOT be given at the beginning of the prefix name.
2980    The first character of all runtime options is AUTOMATICALLY the
2981    hyphen.
2982 
2983    Level: advanced
2984 
2985 .keywords: TS, append, options, prefix, database
2986 
2987 .seealso: TSGetOptionsPrefix()
2988 
2989 @*/
2990 PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2991 {
2992   PetscErrorCode ierr;
2993   SNES           snes;
2994 
2995   PetscFunctionBegin;
2996   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2997   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2998   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2999   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
3000   PetscFunctionReturn(0);
3001 }
3002 
3003 #undef __FUNCT__
3004 #define __FUNCT__ "TSGetOptionsPrefix"
3005 /*@C
3006    TSGetOptionsPrefix - Sets the prefix used for searching for all
3007    TS options in the database.
3008 
3009    Not Collective
3010 
3011    Input Parameter:
3012 .  ts - The TS context
3013 
3014    Output Parameter:
3015 .  prefix - A pointer to the prefix string used
3016 
3017    Notes: On the fortran side, the user should pass in a string 'prifix' of
3018    sufficient length to hold the prefix.
3019 
3020    Level: intermediate
3021 
3022 .keywords: TS, get, options, prefix, database
3023 
3024 .seealso: TSAppendOptionsPrefix()
3025 @*/
3026 PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
3027 {
3028   PetscErrorCode ierr;
3029 
3030   PetscFunctionBegin;
3031   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3032   PetscValidPointer(prefix,2);
3033   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
3034   PetscFunctionReturn(0);
3035 }
3036 
3037 #undef __FUNCT__
3038 #define __FUNCT__ "TSGetRHSJacobian"
3039 /*@C
3040    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
3041 
3042    Not Collective, but parallel objects are returned if TS is parallel
3043 
3044    Input Parameter:
3045 .  ts  - The TS context obtained from TSCreate()
3046 
3047    Output Parameters:
3048 +  Amat - The (approximate) Jacobian J of G, where U_t = G(U,t)  (or NULL)
3049 .  Pmat - The matrix from which the preconditioner is constructed, usually the same as Amat  (or NULL)
3050 .  func - Function to compute the Jacobian of the RHS  (or NULL)
3051 -  ctx - User-defined context for Jacobian evaluation routine  (or NULL)
3052 
3053    Notes: You can pass in NULL for any return argument you do not need.
3054 
3055    Level: intermediate
3056 
3057 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
3058 
3059 .keywords: TS, timestep, get, matrix, Jacobian
3060 @*/
3061 PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *Amat,Mat *Pmat,TSRHSJacobian *func,void **ctx)
3062 {
3063   PetscErrorCode ierr;
3064   SNES           snes;
3065   DM             dm;
3066 
3067   PetscFunctionBegin;
3068   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3069   ierr = SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);CHKERRQ(ierr);
3070   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
3071   ierr = DMTSGetRHSJacobian(dm,func,ctx);CHKERRQ(ierr);
3072   PetscFunctionReturn(0);
3073 }
3074 
3075 #undef __FUNCT__
3076 #define __FUNCT__ "TSGetIJacobian"
3077 /*@C
3078    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
3079 
3080    Not Collective, but parallel objects are returned if TS is parallel
3081 
3082    Input Parameter:
3083 .  ts  - The TS context obtained from TSCreate()
3084 
3085    Output Parameters:
3086 +  Amat  - The (approximate) Jacobian of F(t,U,U_t)
3087 .  Pmat - The matrix from which the preconditioner is constructed, often the same as Amat
3088 .  f   - The function to compute the matrices
3089 - ctx - User-defined context for Jacobian evaluation routine
3090 
3091    Notes: You can pass in NULL for any return argument you do not need.
3092 
3093    Level: advanced
3094 
3095 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
3096 
3097 .keywords: TS, timestep, get, matrix, Jacobian
3098 @*/
3099 PetscErrorCode  TSGetIJacobian(TS ts,Mat *Amat,Mat *Pmat,TSIJacobian *f,void **ctx)
3100 {
3101   PetscErrorCode ierr;
3102   SNES           snes;
3103   DM             dm;
3104 
3105   PetscFunctionBegin;
3106   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3107   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
3108   ierr = SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);CHKERRQ(ierr);
3109   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
3110   ierr = DMTSGetIJacobian(dm,f,ctx);CHKERRQ(ierr);
3111   PetscFunctionReturn(0);
3112 }
3113 
3114 
3115 #undef __FUNCT__
3116 #define __FUNCT__ "TSMonitorDrawSolution"
3117 /*@C
3118    TSMonitorDrawSolution - Monitors progress of the TS solvers by calling
3119    VecView() for the solution at each timestep
3120 
3121    Collective on TS
3122 
3123    Input Parameters:
3124 +  ts - the TS context
3125 .  step - current time-step
3126 .  ptime - current time
3127 -  dummy - either a viewer or NULL
3128 
3129    Options Database:
3130 .   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
3131 
3132    Notes: the initial solution and current solution are not displayed with a common axis scaling so generally the option -ts_monitor_draw_solution_initial
3133        will look bad
3134 
3135    Level: intermediate
3136 
3137 .keywords: TS,  vector, monitor, view
3138 
3139 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3140 @*/
3141 PetscErrorCode  TSMonitorDrawSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
3142 {
3143   PetscErrorCode   ierr;
3144   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
3145   PetscDraw        draw;
3146 
3147   PetscFunctionBegin;
3148   if (!step && ictx->showinitial) {
3149     if (!ictx->initialsolution) {
3150       ierr = VecDuplicate(u,&ictx->initialsolution);CHKERRQ(ierr);
3151     }
3152     ierr = VecCopy(u,ictx->initialsolution);CHKERRQ(ierr);
3153   }
3154   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
3155 
3156   if (ictx->showinitial) {
3157     PetscReal pause;
3158     ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr);
3159     ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr);
3160     ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr);
3161     ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr);
3162     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr);
3163   }
3164   ierr = VecView(u,ictx->viewer);CHKERRQ(ierr);
3165   if (ictx->showtimestepandtime) {
3166     PetscReal xl,yl,xr,yr,tw,w,h;
3167     char      time[32];
3168     size_t    len;
3169 
3170     ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
3171     ierr = PetscSNPrintf(time,32,"Timestep %d Time %f",(int)step,(double)ptime);CHKERRQ(ierr);
3172     ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
3173     ierr =  PetscStrlen(time,&len);CHKERRQ(ierr);
3174     ierr = PetscDrawStringGetSize(draw,&tw,NULL);CHKERRQ(ierr);
3175     w    = xl + .5*(xr - xl) - .5*len*tw;
3176     h    = yl + .95*(yr - yl);
3177     ierr = PetscDrawString(draw,w,h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
3178     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
3179   }
3180 
3181   if (ictx->showinitial) {
3182     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr);
3183   }
3184   PetscFunctionReturn(0);
3185 }
3186 
3187 #undef __FUNCT__
3188 #define __FUNCT__ "TSMonitorDrawSolutionPhase"
3189 /*@C
3190    TSMonitorDrawSolutionPhase - Monitors progress of the TS solvers by plotting the solution as a phase diagram
3191 
3192    Collective on TS
3193 
3194    Input Parameters:
3195 +  ts - the TS context
3196 .  step - current time-step
3197 .  ptime - current time
3198 -  dummy - either a viewer or NULL
3199 
3200    Level: intermediate
3201 
3202 .keywords: TS,  vector, monitor, view
3203 
3204 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3205 @*/
3206 PetscErrorCode  TSMonitorDrawSolutionPhase(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
3207 {
3208   PetscErrorCode    ierr;
3209   TSMonitorDrawCtx  ictx = (TSMonitorDrawCtx)dummy;
3210   PetscDraw         draw;
3211   MPI_Comm          comm;
3212   PetscInt          n;
3213   PetscMPIInt       size;
3214   PetscReal         xl,yl,xr,yr,tw,w,h;
3215   char              time[32];
3216   size_t            len;
3217   const PetscScalar *U;
3218 
3219   PetscFunctionBegin;
3220   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
3221   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
3222   if (size != 1) SETERRQ(comm,PETSC_ERR_SUP,"Only allowed for sequential runs");
3223   ierr = VecGetSize(u,&n);CHKERRQ(ierr);
3224   if (n != 2) SETERRQ(comm,PETSC_ERR_SUP,"Only for ODEs with two unknowns");
3225 
3226   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
3227 
3228   ierr = VecGetArrayRead(u,&U);CHKERRQ(ierr);
3229   ierr = PetscDrawAxisGetLimits(ictx->axis,&xl,&xr,&yl,&yr);CHKERRQ(ierr);
3230   if ((PetscRealPart(U[0]) < xl) || (PetscRealPart(U[1]) < yl) || (PetscRealPart(U[0]) > xr) || (PetscRealPart(U[1]) > yr)) {
3231       ierr = VecRestoreArrayRead(u,&U);CHKERRQ(ierr);
3232       PetscFunctionReturn(0);
3233   }
3234   if (!step) ictx->color++;
3235   ierr = PetscDrawPoint(draw,PetscRealPart(U[0]),PetscRealPart(U[1]),ictx->color);CHKERRQ(ierr);
3236   ierr = VecRestoreArrayRead(u,&U);CHKERRQ(ierr);
3237 
3238   if (ictx->showtimestepandtime) {
3239     ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
3240     ierr = PetscSNPrintf(time,32,"Timestep %d Time %f",(int)step,(double)ptime);CHKERRQ(ierr);
3241     ierr = PetscStrlen(time,&len);CHKERRQ(ierr);
3242     ierr = PetscDrawStringGetSize(draw,&tw,NULL);CHKERRQ(ierr);
3243     w    = xl + .5*(xr - xl) - .5*len*tw;
3244     h    = yl + .95*(yr - yl);
3245     ierr = PetscDrawString(draw,w,h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
3246   }
3247   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
3248   PetscFunctionReturn(0);
3249 }
3250 
3251 
3252 #undef __FUNCT__
3253 #define __FUNCT__ "TSMonitorDrawCtxDestroy"
3254 /*@C
3255    TSMonitorDrawCtxDestroy - Destroys the monitor context for TSMonitorDrawSolution()
3256 
3257    Collective on TS
3258 
3259    Input Parameters:
3260 .    ctx - the monitor context
3261 
3262    Level: intermediate
3263 
3264 .keywords: TS,  vector, monitor, view
3265 
3266 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawSolution(), TSMonitorDrawError()
3267 @*/
3268 PetscErrorCode  TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
3269 {
3270   PetscErrorCode ierr;
3271 
3272   PetscFunctionBegin;
3273   ierr = PetscDrawAxisDestroy(&(*ictx)->axis);CHKERRQ(ierr);
3274   ierr = PetscViewerDestroy(&(*ictx)->viewer);CHKERRQ(ierr);
3275   ierr = VecDestroy(&(*ictx)->initialsolution);CHKERRQ(ierr);
3276   ierr = PetscFree(*ictx);CHKERRQ(ierr);
3277   PetscFunctionReturn(0);
3278 }
3279 
3280 #undef __FUNCT__
3281 #define __FUNCT__ "TSMonitorDrawCtxCreate"
3282 /*@C
3283    TSMonitorDrawCtxCreate - Creates the monitor context for TSMonitorDrawCtx
3284 
3285    Collective on TS
3286 
3287    Input Parameter:
3288 .    ts - time-step context
3289 
3290    Output Patameter:
3291 .    ctx - the monitor context
3292 
3293    Options Database:
3294 .   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
3295 
3296    Level: intermediate
3297 
3298 .keywords: TS,  vector, monitor, view
3299 
3300 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawCtx()
3301 @*/
3302 PetscErrorCode  TSMonitorDrawCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorDrawCtx *ctx)
3303 {
3304   PetscErrorCode   ierr;
3305 
3306   PetscFunctionBegin;
3307   ierr = PetscNew(struct _n_TSMonitorDrawCtx,ctx);CHKERRQ(ierr);
3308   ierr = PetscViewerDrawOpen(comm,host,label,x,y,m,n,&(*ctx)->viewer);CHKERRQ(ierr);
3309   ierr = PetscViewerSetFromOptions((*ctx)->viewer);CHKERRQ(ierr);
3310 
3311   (*ctx)->howoften    = howoften;
3312   (*ctx)->showinitial = PETSC_FALSE;
3313   ierr = PetscOptionsGetBool(NULL,"-ts_monitor_draw_solution_initial",&(*ctx)->showinitial,NULL);CHKERRQ(ierr);
3314 
3315   (*ctx)->showtimestepandtime = PETSC_FALSE;
3316   ierr = PetscOptionsGetBool(NULL,"-ts_monitor_draw_solution_show_time",&(*ctx)->showtimestepandtime,NULL);CHKERRQ(ierr);
3317   (*ctx)->color = PETSC_DRAW_WHITE;
3318   PetscFunctionReturn(0);
3319 }
3320 
3321 #undef __FUNCT__
3322 #define __FUNCT__ "TSMonitorDrawError"
3323 /*@C
3324    TSMonitorDrawError - Monitors progress of the TS solvers by calling
3325    VecView() for the error at each timestep
3326 
3327    Collective on TS
3328 
3329    Input Parameters:
3330 +  ts - the TS context
3331 .  step - current time-step
3332 .  ptime - current time
3333 -  dummy - either a viewer or NULL
3334 
3335    Level: intermediate
3336 
3337 .keywords: TS,  vector, monitor, view
3338 
3339 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3340 @*/
3341 PetscErrorCode  TSMonitorDrawError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
3342 {
3343   PetscErrorCode   ierr;
3344   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
3345   PetscViewer      viewer = ctx->viewer;
3346   Vec              work;
3347 
3348   PetscFunctionBegin;
3349   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
3350   ierr = VecDuplicate(u,&work);CHKERRQ(ierr);
3351   ierr = TSComputeSolutionFunction(ts,ptime,work);CHKERRQ(ierr);
3352   ierr = VecAXPY(work,-1.0,u);CHKERRQ(ierr);
3353   ierr = VecView(work,viewer);CHKERRQ(ierr);
3354   ierr = VecDestroy(&work);CHKERRQ(ierr);
3355   PetscFunctionReturn(0);
3356 }
3357 
3358 #include <petsc-private/dmimpl.h>
3359 #undef __FUNCT__
3360 #define __FUNCT__ "TSSetDM"
3361 /*@
3362    TSSetDM - Sets the DM that may be used by some preconditioners
3363 
3364    Logically Collective on TS and DM
3365 
3366    Input Parameters:
3367 +  ts - the preconditioner context
3368 -  dm - the dm
3369 
3370    Level: intermediate
3371 
3372 
3373 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
3374 @*/
3375 PetscErrorCode  TSSetDM(TS ts,DM dm)
3376 {
3377   PetscErrorCode ierr;
3378   SNES           snes;
3379   DMTS           tsdm;
3380 
3381   PetscFunctionBegin;
3382   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3383   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
3384   if (ts->dm) {               /* Move the DMTS context over to the new DM unless the new DM already has one */
3385     if (ts->dm->dmts && !dm->dmts) {
3386       ierr = DMCopyDMTS(ts->dm,dm);CHKERRQ(ierr);
3387       ierr = DMGetDMTS(ts->dm,&tsdm);CHKERRQ(ierr);
3388       if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
3389         tsdm->originaldm = dm;
3390       }
3391     }
3392     ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
3393   }
3394   ts->dm = dm;
3395 
3396   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3397   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
3398   PetscFunctionReturn(0);
3399 }
3400 
3401 #undef __FUNCT__
3402 #define __FUNCT__ "TSGetDM"
3403 /*@
3404    TSGetDM - Gets the DM that may be used by some preconditioners
3405 
3406    Not Collective
3407 
3408    Input Parameter:
3409 . ts - the preconditioner context
3410 
3411    Output Parameter:
3412 .  dm - the dm
3413 
3414    Level: intermediate
3415 
3416 
3417 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
3418 @*/
3419 PetscErrorCode  TSGetDM(TS ts,DM *dm)
3420 {
3421   PetscErrorCode ierr;
3422 
3423   PetscFunctionBegin;
3424   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3425   if (!ts->dm) {
3426     ierr = DMShellCreate(PetscObjectComm((PetscObject)ts),&ts->dm);CHKERRQ(ierr);
3427     if (ts->snes) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
3428   }
3429   *dm = ts->dm;
3430   PetscFunctionReturn(0);
3431 }
3432 
3433 #undef __FUNCT__
3434 #define __FUNCT__ "SNESTSFormFunction"
3435 /*@
3436    SNESTSFormFunction - Function to evaluate nonlinear residual
3437 
3438    Logically Collective on SNES
3439 
3440    Input Parameter:
3441 + snes - nonlinear solver
3442 . U - the current state at which to evaluate the residual
3443 - ctx - user context, must be a TS
3444 
3445    Output Parameter:
3446 . F - the nonlinear residual
3447 
3448    Notes:
3449    This function is not normally called by users and is automatically registered with the SNES used by TS.
3450    It is most frequently passed to MatFDColoringSetFunction().
3451 
3452    Level: advanced
3453 
3454 .seealso: SNESSetFunction(), MatFDColoringSetFunction()
3455 @*/
3456 PetscErrorCode  SNESTSFormFunction(SNES snes,Vec U,Vec F,void *ctx)
3457 {
3458   TS             ts = (TS)ctx;
3459   PetscErrorCode ierr;
3460 
3461   PetscFunctionBegin;
3462   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3463   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
3464   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
3465   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
3466   ierr = (ts->ops->snesfunction)(snes,U,F,ts);CHKERRQ(ierr);
3467   PetscFunctionReturn(0);
3468 }
3469 
3470 #undef __FUNCT__
3471 #define __FUNCT__ "SNESTSFormJacobian"
3472 /*@
3473    SNESTSFormJacobian - Function to evaluate the Jacobian
3474 
3475    Collective on SNES
3476 
3477    Input Parameter:
3478 + snes - nonlinear solver
3479 . U - the current state at which to evaluate the residual
3480 - ctx - user context, must be a TS
3481 
3482    Output Parameter:
3483 + A - the Jacobian
3484 . B - the preconditioning matrix (may be the same as A)
3485 - flag - indicates any structure change in the matrix
3486 
3487    Notes:
3488    This function is not normally called by users and is automatically registered with the SNES used by TS.
3489 
3490    Level: developer
3491 
3492 .seealso: SNESSetJacobian()
3493 @*/
3494 PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec U,Mat *A,Mat *B,MatStructure *flag,void *ctx)
3495 {
3496   TS             ts = (TS)ctx;
3497   PetscErrorCode ierr;
3498 
3499   PetscFunctionBegin;
3500   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3501   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
3502   PetscValidPointer(A,3);
3503   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
3504   PetscValidPointer(B,4);
3505   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
3506   PetscValidPointer(flag,5);
3507   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
3508   ierr = (ts->ops->snesjacobian)(snes,U,A,B,flag,ts);CHKERRQ(ierr);
3509   PetscFunctionReturn(0);
3510 }
3511 
3512 #undef __FUNCT__
3513 #define __FUNCT__ "TSComputeRHSFunctionLinear"
3514 /*@C
3515    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only
3516 
3517    Collective on TS
3518 
3519    Input Arguments:
3520 +  ts - time stepping context
3521 .  t - time at which to evaluate
3522 .  U - state at which to evaluate
3523 -  ctx - context
3524 
3525    Output Arguments:
3526 .  F - right hand side
3527 
3528    Level: intermediate
3529 
3530    Notes:
3531    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
3532    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
3533 
3534 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
3535 @*/
3536 PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
3537 {
3538   PetscErrorCode ierr;
3539   Mat            Arhs,Brhs;
3540   MatStructure   flg2;
3541 
3542   PetscFunctionBegin;
3543   ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
3544   ierr = TSComputeRHSJacobian(ts,t,U,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
3545   ierr = MatMult(Arhs,U,F);CHKERRQ(ierr);
3546   PetscFunctionReturn(0);
3547 }
3548 
3549 #undef __FUNCT__
3550 #define __FUNCT__ "TSComputeRHSJacobianConstant"
3551 /*@C
3552    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
3553 
3554    Collective on TS
3555 
3556    Input Arguments:
3557 +  ts - time stepping context
3558 .  t - time at which to evaluate
3559 .  U - state at which to evaluate
3560 -  ctx - context
3561 
3562    Output Arguments:
3563 +  A - pointer to operator
3564 .  B - pointer to preconditioning matrix
3565 -  flg - matrix structure flag
3566 
3567    Level: intermediate
3568 
3569    Notes:
3570    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
3571 
3572 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
3573 @*/
3574 PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec U,Mat *A,Mat *B,MatStructure *flg,void *ctx)
3575 {
3576   PetscFunctionBegin;
3577   *flg = SAME_PRECONDITIONER;
3578   PetscFunctionReturn(0);
3579 }
3580 
3581 #undef __FUNCT__
3582 #define __FUNCT__ "TSComputeIFunctionLinear"
3583 /*@C
3584    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
3585 
3586    Collective on TS
3587 
3588    Input Arguments:
3589 +  ts - time stepping context
3590 .  t - time at which to evaluate
3591 .  U - state at which to evaluate
3592 .  Udot - time derivative of state vector
3593 -  ctx - context
3594 
3595    Output Arguments:
3596 .  F - left hand side
3597 
3598    Level: intermediate
3599 
3600    Notes:
3601    The assumption here is that the left hand side is of the form A*Udot (and not A*Udot + B*U). For other cases, the
3602    user is required to write their own TSComputeIFunction.
3603    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
3604    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
3605 
3606 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
3607 @*/
3608 PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
3609 {
3610   PetscErrorCode ierr;
3611   Mat            A,B;
3612   MatStructure   flg2;
3613 
3614   PetscFunctionBegin;
3615   ierr = TSGetIJacobian(ts,&A,&B,NULL,NULL);CHKERRQ(ierr);
3616   ierr = TSComputeIJacobian(ts,t,U,Udot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr);
3617   ierr = MatMult(A,Udot,F);CHKERRQ(ierr);
3618   PetscFunctionReturn(0);
3619 }
3620 
3621 #undef __FUNCT__
3622 #define __FUNCT__ "TSComputeIJacobianConstant"
3623 /*@C
3624    TSComputeIJacobianConstant - Reuses a time-independent for a semi-implicit DAE or ODE
3625 
3626    Collective on TS
3627 
3628    Input Arguments:
3629 +  ts - time stepping context
3630 .  t - time at which to evaluate
3631 .  U - state at which to evaluate
3632 .  Udot - time derivative of state vector
3633 .  shift - shift to apply
3634 -  ctx - context
3635 
3636    Output Arguments:
3637 +  A - pointer to operator
3638 .  B - pointer to preconditioning matrix
3639 -  flg - matrix structure flag
3640 
3641    Level: advanced
3642 
3643    Notes:
3644    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
3645 
3646    It is only appropriate for problems of the form
3647 
3648 $     M Udot = F(U,t)
3649 
3650   where M is constant and F is non-stiff.  The user must pass M to TSSetIJacobian().  The current implementation only
3651   works with IMEX time integration methods such as TSROSW and TSARKIMEX, since there is no support for de-constructing
3652   an implicit operator of the form
3653 
3654 $    shift*M + J
3655 
3656   where J is the Jacobian of -F(U).  Support may be added in a future version of PETSc, but for now, the user must store
3657   a copy of M or reassemble it when requested.
3658 
3659 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
3660 @*/
3661 PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
3662 {
3663   PetscErrorCode ierr;
3664 
3665   PetscFunctionBegin;
3666   ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr);
3667   ts->ijacobian.shift = shift;
3668   *flg = SAME_PRECONDITIONER;
3669   PetscFunctionReturn(0);
3670 }
3671 
3672 #undef __FUNCT__
3673 #define __FUNCT__ "TSGetEquationType"
3674 /*@
3675    TSGetEquationType - Gets the type of the equation that TS is solving.
3676 
3677    Not Collective
3678 
3679    Input Parameter:
3680 .  ts - the TS context
3681 
3682    Output Parameter:
3683 .  equation_type - see TSEquationType
3684 
3685    Level: beginner
3686 
3687 .keywords: TS, equation type
3688 
3689 .seealso: TSSetEquationType(), TSEquationType
3690 @*/
3691 PetscErrorCode  TSGetEquationType(TS ts,TSEquationType *equation_type)
3692 {
3693   PetscFunctionBegin;
3694   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3695   PetscValidPointer(equation_type,2);
3696   *equation_type = ts->equation_type;
3697   PetscFunctionReturn(0);
3698 }
3699 
3700 #undef __FUNCT__
3701 #define __FUNCT__ "TSSetEquationType"
3702 /*@
3703    TSSetEquationType - Sets the type of the equation that TS is solving.
3704 
3705    Not Collective
3706 
3707    Input Parameter:
3708 +  ts - the TS context
3709 .  equation_type - see TSEquationType
3710 
3711    Level: advanced
3712 
3713 .keywords: TS, equation type
3714 
3715 .seealso: TSGetEquationType(), TSEquationType
3716 @*/
3717 PetscErrorCode  TSSetEquationType(TS ts,TSEquationType equation_type)
3718 {
3719   PetscFunctionBegin;
3720   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3721   ts->equation_type = equation_type;
3722   PetscFunctionReturn(0);
3723 }
3724 
3725 #undef __FUNCT__
3726 #define __FUNCT__ "TSGetConvergedReason"
3727 /*@
3728    TSGetConvergedReason - Gets the reason the TS iteration was stopped.
3729 
3730    Not Collective
3731 
3732    Input Parameter:
3733 .  ts - the TS context
3734 
3735    Output Parameter:
3736 .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
3737             manual pages for the individual convergence tests for complete lists
3738 
3739    Level: beginner
3740 
3741    Notes:
3742    Can only be called after the call to TSSolve() is complete.
3743 
3744 .keywords: TS, nonlinear, set, convergence, test
3745 
3746 .seealso: TSSetConvergenceTest(), TSConvergedReason
3747 @*/
3748 PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
3749 {
3750   PetscFunctionBegin;
3751   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3752   PetscValidPointer(reason,2);
3753   *reason = ts->reason;
3754   PetscFunctionReturn(0);
3755 }
3756 
3757 #undef __FUNCT__
3758 #define __FUNCT__ "TSSetConvergedReason"
3759 /*@
3760    TSSetConvergedReason - Sets the reason for handling the convergence of TSSolve.
3761 
3762    Not Collective
3763 
3764    Input Parameter:
3765 +  ts - the TS context
3766 .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
3767             manual pages for the individual convergence tests for complete lists
3768 
3769    Level: advanced
3770 
3771    Notes:
3772    Can only be called during TSSolve() is active.
3773 
3774 .keywords: TS, nonlinear, set, convergence, test
3775 
3776 .seealso: TSConvergedReason
3777 @*/
3778 PetscErrorCode  TSSetConvergedReason(TS ts,TSConvergedReason reason)
3779 {
3780   PetscFunctionBegin;
3781   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3782   ts->reason = reason;
3783   PetscFunctionReturn(0);
3784 }
3785 
3786 #undef __FUNCT__
3787 #define __FUNCT__ "TSGetSolveTime"
3788 /*@
3789    TSGetSolveTime - Gets the time after a call to TSSolve()
3790 
3791    Not Collective
3792 
3793    Input Parameter:
3794 .  ts - the TS context
3795 
3796    Output Parameter:
3797 .  ftime - the final time. This time should correspond to the final time set with TSSetDuration()
3798 
3799    Level: beginner
3800 
3801    Notes:
3802    Can only be called after the call to TSSolve() is complete.
3803 
3804 .keywords: TS, nonlinear, set, convergence, test
3805 
3806 .seealso: TSSetConvergenceTest(), TSConvergedReason
3807 @*/
3808 PetscErrorCode  TSGetSolveTime(TS ts,PetscReal *ftime)
3809 {
3810   PetscFunctionBegin;
3811   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3812   PetscValidPointer(ftime,2);
3813   *ftime = ts->solvetime;
3814   PetscFunctionReturn(0);
3815 }
3816 
3817 #undef __FUNCT__
3818 #define __FUNCT__ "TSGetSNESIterations"
3819 /*@
3820    TSGetSNESIterations - Gets the total number of nonlinear iterations
3821    used by the time integrator.
3822 
3823    Not Collective
3824 
3825    Input Parameter:
3826 .  ts - TS context
3827 
3828    Output Parameter:
3829 .  nits - number of nonlinear iterations
3830 
3831    Notes:
3832    This counter is reset to zero for each successive call to TSSolve().
3833 
3834    Level: intermediate
3835 
3836 .keywords: TS, get, number, nonlinear, iterations
3837 
3838 .seealso:  TSGetKSPIterations()
3839 @*/
3840 PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
3841 {
3842   PetscFunctionBegin;
3843   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3844   PetscValidIntPointer(nits,2);
3845   *nits = ts->snes_its;
3846   PetscFunctionReturn(0);
3847 }
3848 
3849 #undef __FUNCT__
3850 #define __FUNCT__ "TSGetKSPIterations"
3851 /*@
3852    TSGetKSPIterations - Gets the total number of linear iterations
3853    used by the time integrator.
3854 
3855    Not Collective
3856 
3857    Input Parameter:
3858 .  ts - TS context
3859 
3860    Output Parameter:
3861 .  lits - number of linear iterations
3862 
3863    Notes:
3864    This counter is reset to zero for each successive call to TSSolve().
3865 
3866    Level: intermediate
3867 
3868 .keywords: TS, get, number, linear, iterations
3869 
3870 .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
3871 @*/
3872 PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
3873 {
3874   PetscFunctionBegin;
3875   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3876   PetscValidIntPointer(lits,2);
3877   *lits = ts->ksp_its;
3878   PetscFunctionReturn(0);
3879 }
3880 
3881 #undef __FUNCT__
3882 #define __FUNCT__ "TSGetStepRejections"
3883 /*@
3884    TSGetStepRejections - Gets the total number of rejected steps.
3885 
3886    Not Collective
3887 
3888    Input Parameter:
3889 .  ts - TS context
3890 
3891    Output Parameter:
3892 .  rejects - number of steps rejected
3893 
3894    Notes:
3895    This counter is reset to zero for each successive call to TSSolve().
3896 
3897    Level: intermediate
3898 
3899 .keywords: TS, get, number
3900 
3901 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
3902 @*/
3903 PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
3904 {
3905   PetscFunctionBegin;
3906   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3907   PetscValidIntPointer(rejects,2);
3908   *rejects = ts->reject;
3909   PetscFunctionReturn(0);
3910 }
3911 
3912 #undef __FUNCT__
3913 #define __FUNCT__ "TSGetSNESFailures"
3914 /*@
3915    TSGetSNESFailures - Gets the total number of failed SNES solves
3916 
3917    Not Collective
3918 
3919    Input Parameter:
3920 .  ts - TS context
3921 
3922    Output Parameter:
3923 .  fails - number of failed nonlinear solves
3924 
3925    Notes:
3926    This counter is reset to zero for each successive call to TSSolve().
3927 
3928    Level: intermediate
3929 
3930 .keywords: TS, get, number
3931 
3932 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
3933 @*/
3934 PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
3935 {
3936   PetscFunctionBegin;
3937   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3938   PetscValidIntPointer(fails,2);
3939   *fails = ts->num_snes_failures;
3940   PetscFunctionReturn(0);
3941 }
3942 
3943 #undef __FUNCT__
3944 #define __FUNCT__ "TSSetMaxStepRejections"
3945 /*@
3946    TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails
3947 
3948    Not Collective
3949 
3950    Input Parameter:
3951 +  ts - TS context
3952 -  rejects - maximum number of rejected steps, pass -1 for unlimited
3953 
3954    Notes:
3955    The counter is reset to zero for each step
3956 
3957    Options Database Key:
3958  .  -ts_max_reject - Maximum number of step rejections before a step fails
3959 
3960    Level: intermediate
3961 
3962 .keywords: TS, set, maximum, number
3963 
3964 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3965 @*/
3966 PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
3967 {
3968   PetscFunctionBegin;
3969   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3970   ts->max_reject = rejects;
3971   PetscFunctionReturn(0);
3972 }
3973 
3974 #undef __FUNCT__
3975 #define __FUNCT__ "TSSetMaxSNESFailures"
3976 /*@
3977    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves
3978 
3979    Not Collective
3980 
3981    Input Parameter:
3982 +  ts - TS context
3983 -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited
3984 
3985    Notes:
3986    The counter is reset to zero for each successive call to TSSolve().
3987 
3988    Options Database Key:
3989  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures
3990 
3991    Level: intermediate
3992 
3993 .keywords: TS, set, maximum, number
3994 
3995 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
3996 @*/
3997 PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
3998 {
3999   PetscFunctionBegin;
4000   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4001   ts->max_snes_failures = fails;
4002   PetscFunctionReturn(0);
4003 }
4004 
4005 #undef __FUNCT__
4006 #define __FUNCT__ "TSSetErrorIfStepFails()"
4007 /*@
4008    TSSetErrorIfStepFails - Error if no step succeeds
4009 
4010    Not Collective
4011 
4012    Input Parameter:
4013 +  ts - TS context
4014 -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure
4015 
4016    Options Database Key:
4017  .  -ts_error_if_step_fails - Error if no step succeeds
4018 
4019    Level: intermediate
4020 
4021 .keywords: TS, set, error
4022 
4023 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
4024 @*/
4025 PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
4026 {
4027   PetscFunctionBegin;
4028   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4029   ts->errorifstepfailed = err;
4030   PetscFunctionReturn(0);
4031 }
4032 
4033 #undef __FUNCT__
4034 #define __FUNCT__ "TSMonitorSolutionBinary"
4035 /*@C
4036    TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file
4037 
4038    Collective on TS
4039 
4040    Input Parameters:
4041 +  ts - the TS context
4042 .  step - current time-step
4043 .  ptime - current time
4044 .  u - current state
4045 -  viewer - binary viewer
4046 
4047    Level: intermediate
4048 
4049 .keywords: TS,  vector, monitor, view
4050 
4051 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4052 @*/
4053 PetscErrorCode  TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec u,void *viewer)
4054 {
4055   PetscErrorCode ierr;
4056   PetscViewer    v = (PetscViewer)viewer;
4057 
4058   PetscFunctionBegin;
4059   ierr = VecView(u,v);CHKERRQ(ierr);
4060   PetscFunctionReturn(0);
4061 }
4062 
4063 #undef __FUNCT__
4064 #define __FUNCT__ "TSMonitorSolutionVTK"
4065 /*@C
4066    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.
4067 
4068    Collective on TS
4069 
4070    Input Parameters:
4071 +  ts - the TS context
4072 .  step - current time-step
4073 .  ptime - current time
4074 .  u - current state
4075 -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
4076 
4077    Level: intermediate
4078 
4079    Notes:
4080    The VTK format does not allow writing multiple time steps in the same file, therefore a different file will be written for each time step.
4081    These are named according to the file name template.
4082 
4083    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().
4084 
4085 .keywords: TS,  vector, monitor, view
4086 
4087 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4088 @*/
4089 PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec u,void *filenametemplate)
4090 {
4091   PetscErrorCode ierr;
4092   char           filename[PETSC_MAX_PATH_LEN];
4093   PetscViewer    viewer;
4094 
4095   PetscFunctionBegin;
4096   ierr = PetscSNPrintf(filename,sizeof(filename),(const char*)filenametemplate,step);CHKERRQ(ierr);
4097   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts),filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
4098   ierr = VecView(u,viewer);CHKERRQ(ierr);
4099   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
4100   PetscFunctionReturn(0);
4101 }
4102 
4103 #undef __FUNCT__
4104 #define __FUNCT__ "TSMonitorSolutionVTKDestroy"
4105 /*@C
4106    TSMonitorSolutionVTKDestroy - Destroy context for monitoring
4107 
4108    Collective on TS
4109 
4110    Input Parameters:
4111 .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
4112 
4113    Level: intermediate
4114 
4115    Note:
4116    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().
4117 
4118 .keywords: TS,  vector, monitor, view
4119 
4120 .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
4121 @*/
4122 PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
4123 {
4124   PetscErrorCode ierr;
4125 
4126   PetscFunctionBegin;
4127   ierr = PetscFree(*(char**)filenametemplate);CHKERRQ(ierr);
4128   PetscFunctionReturn(0);
4129 }
4130 
4131 #undef __FUNCT__
4132 #define __FUNCT__ "TSGetAdapt"
4133 /*@
4134    TSGetAdapt - Get the adaptive controller context for the current method
4135 
4136    Collective on TS if controller has not been created yet
4137 
4138    Input Arguments:
4139 .  ts - time stepping context
4140 
4141    Output Arguments:
4142 .  adapt - adaptive controller
4143 
4144    Level: intermediate
4145 
4146 .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
4147 @*/
4148 PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
4149 {
4150   PetscErrorCode ierr;
4151 
4152   PetscFunctionBegin;
4153   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4154   PetscValidPointer(adapt,2);
4155   if (!ts->adapt) {
4156     ierr = TSAdaptCreate(PetscObjectComm((PetscObject)ts),&ts->adapt);CHKERRQ(ierr);
4157     ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->adapt);CHKERRQ(ierr);
4158     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);CHKERRQ(ierr);
4159   }
4160   *adapt = ts->adapt;
4161   PetscFunctionReturn(0);
4162 }
4163 
4164 #undef __FUNCT__
4165 #define __FUNCT__ "TSSetTolerances"
4166 /*@
4167    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller
4168 
4169    Logically Collective
4170 
4171    Input Arguments:
4172 +  ts - time integration context
4173 .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
4174 .  vatol - vector of absolute tolerances or NULL, used in preference to atol if present
4175 .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
4176 -  vrtol - vector of relative tolerances or NULL, used in preference to atol if present
4177 
4178    Level: beginner
4179 
4180 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
4181 @*/
4182 PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
4183 {
4184   PetscErrorCode ierr;
4185 
4186   PetscFunctionBegin;
4187   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
4188   if (vatol) {
4189     ierr = PetscObjectReference((PetscObject)vatol);CHKERRQ(ierr);
4190     ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
4191 
4192     ts->vatol = vatol;
4193   }
4194   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
4195   if (vrtol) {
4196     ierr = PetscObjectReference((PetscObject)vrtol);CHKERRQ(ierr);
4197     ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
4198 
4199     ts->vrtol = vrtol;
4200   }
4201   PetscFunctionReturn(0);
4202 }
4203 
4204 #undef __FUNCT__
4205 #define __FUNCT__ "TSGetTolerances"
4206 /*@
4207    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
4208 
4209    Logically Collective
4210 
4211    Input Arguments:
4212 .  ts - time integration context
4213 
4214    Output Arguments:
4215 +  atol - scalar absolute tolerances, NULL to ignore
4216 .  vatol - vector of absolute tolerances, NULL to ignore
4217 .  rtol - scalar relative tolerances, NULL to ignore
4218 -  vrtol - vector of relative tolerances, NULL to ignore
4219 
4220    Level: beginner
4221 
4222 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
4223 @*/
4224 PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
4225 {
4226   PetscFunctionBegin;
4227   if (atol)  *atol  = ts->atol;
4228   if (vatol) *vatol = ts->vatol;
4229   if (rtol)  *rtol  = ts->rtol;
4230   if (vrtol) *vrtol = ts->vrtol;
4231   PetscFunctionReturn(0);
4232 }
4233 
4234 #undef __FUNCT__
4235 #define __FUNCT__ "TSErrorNormWRMS"
4236 /*@
4237    TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state
4238 
4239    Collective on TS
4240 
4241    Input Arguments:
4242 +  ts - time stepping context
4243 -  Y - state vector to be compared to ts->vec_sol
4244 
4245    Output Arguments:
4246 .  norm - weighted norm, a value of 1.0 is considered small
4247 
4248    Level: developer
4249 
4250 .seealso: TSSetTolerances()
4251 @*/
4252 PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
4253 {
4254   PetscErrorCode    ierr;
4255   PetscInt          i,n,N;
4256   const PetscScalar *u,*y;
4257   Vec               U;
4258   PetscReal         sum,gsum;
4259 
4260   PetscFunctionBegin;
4261   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4262   PetscValidHeaderSpecific(Y,VEC_CLASSID,2);
4263   PetscValidPointer(norm,3);
4264   U = ts->vec_sol;
4265   PetscCheckSameTypeAndComm(U,1,Y,2);
4266   if (U == Y) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector");
4267 
4268   ierr = VecGetSize(U,&N);CHKERRQ(ierr);
4269   ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr);
4270   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
4271   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
4272   sum  = 0.;
4273   if (ts->vatol && ts->vrtol) {
4274     const PetscScalar *atol,*rtol;
4275     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
4276     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
4277     for (i=0; i<n; i++) {
4278       PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4279       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
4280     }
4281     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
4282     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
4283   } else if (ts->vatol) {       /* vector atol, scalar rtol */
4284     const PetscScalar *atol;
4285     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
4286     for (i=0; i<n; i++) {
4287       PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4288       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
4289     }
4290     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
4291   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
4292     const PetscScalar *rtol;
4293     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
4294     for (i=0; i<n; i++) {
4295       PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4296       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
4297     }
4298     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
4299   } else {                      /* scalar atol, scalar rtol */
4300     for (i=0; i<n; i++) {
4301       PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
4302       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
4303     }
4304   }
4305   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
4306   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
4307 
4308   ierr  = MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
4309   *norm = PetscSqrtReal(gsum / N);
4310   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
4311   PetscFunctionReturn(0);
4312 }
4313 
4314 #undef __FUNCT__
4315 #define __FUNCT__ "TSSetCFLTimeLocal"
4316 /*@
4317    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
4318 
4319    Logically Collective on TS
4320 
4321    Input Arguments:
4322 +  ts - time stepping context
4323 -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)
4324 
4325    Note:
4326    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
4327 
4328    Level: intermediate
4329 
4330 .seealso: TSGetCFLTime(), TSADAPTCFL
4331 @*/
4332 PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
4333 {
4334   PetscFunctionBegin;
4335   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4336   ts->cfltime_local = cfltime;
4337   ts->cfltime       = -1.;
4338   PetscFunctionReturn(0);
4339 }
4340 
4341 #undef __FUNCT__
4342 #define __FUNCT__ "TSGetCFLTime"
4343 /*@
4344    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
4345 
4346    Collective on TS
4347 
4348    Input Arguments:
4349 .  ts - time stepping context
4350 
4351    Output Arguments:
4352 .  cfltime - maximum stable time step for forward Euler
4353 
4354    Level: advanced
4355 
4356 .seealso: TSSetCFLTimeLocal()
4357 @*/
4358 PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
4359 {
4360   PetscErrorCode ierr;
4361 
4362   PetscFunctionBegin;
4363   if (ts->cfltime < 0) {
4364     ierr = MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
4365   }
4366   *cfltime = ts->cfltime;
4367   PetscFunctionReturn(0);
4368 }
4369 
4370 #undef __FUNCT__
4371 #define __FUNCT__ "TSVISetVariableBounds"
4372 /*@
4373    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
4374 
4375    Input Parameters:
4376 .  ts   - the TS context.
4377 .  xl   - lower bound.
4378 .  xu   - upper bound.
4379 
4380    Notes:
4381    If this routine is not called then the lower and upper bounds are set to
4382    SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp().
4383 
4384    Level: advanced
4385 
4386 @*/
4387 PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
4388 {
4389   PetscErrorCode ierr;
4390   SNES           snes;
4391 
4392   PetscFunctionBegin;
4393   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
4394   ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr);
4395   PetscFunctionReturn(0);
4396 }
4397 
4398 #if defined(PETSC_HAVE_MATLAB_ENGINE)
4399 #include <mex.h>
4400 
4401 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
4402 
4403 #undef __FUNCT__
4404 #define __FUNCT__ "TSComputeFunction_Matlab"
4405 /*
4406    TSComputeFunction_Matlab - Calls the function that has been set with
4407                          TSSetFunctionMatlab().
4408 
4409    Collective on TS
4410 
4411    Input Parameters:
4412 +  snes - the TS context
4413 -  u - input vector
4414 
4415    Output Parameter:
4416 .  y - function vector, as set by TSSetFunction()
4417 
4418    Notes:
4419    TSComputeFunction() is typically used within nonlinear solvers
4420    implementations, so most users would not generally call this routine
4421    themselves.
4422 
4423    Level: developer
4424 
4425 .keywords: TS, nonlinear, compute, function
4426 
4427 .seealso: TSSetFunction(), TSGetFunction()
4428 */
4429 PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec u,Vec udot,Vec y, void *ctx)
4430 {
4431   PetscErrorCode  ierr;
4432   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
4433   int             nlhs  = 1,nrhs = 7;
4434   mxArray         *plhs[1],*prhs[7];
4435   long long int   lx = 0,lxdot = 0,ly = 0,ls = 0;
4436 
4437   PetscFunctionBegin;
4438   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
4439   PetscValidHeaderSpecific(u,VEC_CLASSID,3);
4440   PetscValidHeaderSpecific(udot,VEC_CLASSID,4);
4441   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
4442   PetscCheckSameComm(snes,1,u,3);
4443   PetscCheckSameComm(snes,1,y,5);
4444 
4445   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
4446   ierr = PetscMemcpy(&lx,&u,sizeof(u));CHKERRQ(ierr);
4447   ierr = PetscMemcpy(&lxdot,&udot,sizeof(udot));CHKERRQ(ierr);
4448   ierr = PetscMemcpy(&ly,&y,sizeof(u));CHKERRQ(ierr);
4449 
4450   prhs[0] =  mxCreateDoubleScalar((double)ls);
4451   prhs[1] =  mxCreateDoubleScalar(time);
4452   prhs[2] =  mxCreateDoubleScalar((double)lx);
4453   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
4454   prhs[4] =  mxCreateDoubleScalar((double)ly);
4455   prhs[5] =  mxCreateString(sctx->funcname);
4456   prhs[6] =  sctx->ctx;
4457   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
4458   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4459   mxDestroyArray(prhs[0]);
4460   mxDestroyArray(prhs[1]);
4461   mxDestroyArray(prhs[2]);
4462   mxDestroyArray(prhs[3]);
4463   mxDestroyArray(prhs[4]);
4464   mxDestroyArray(prhs[5]);
4465   mxDestroyArray(plhs[0]);
4466   PetscFunctionReturn(0);
4467 }
4468 
4469 
4470 #undef __FUNCT__
4471 #define __FUNCT__ "TSSetFunctionMatlab"
4472 /*
4473    TSSetFunctionMatlab - Sets the function evaluation routine and function
4474    vector for use by the TS routines in solving ODEs
4475    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
4476 
4477    Logically Collective on TS
4478 
4479    Input Parameters:
4480 +  ts - the TS context
4481 -  func - function evaluation routine
4482 
4483    Calling sequence of func:
4484 $    func (TS ts,PetscReal time,Vec u,Vec udot,Vec f,void *ctx);
4485 
4486    Level: beginner
4487 
4488 .keywords: TS, nonlinear, set, function
4489 
4490 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
4491 */
4492 PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
4493 {
4494   PetscErrorCode  ierr;
4495   TSMatlabContext *sctx;
4496 
4497   PetscFunctionBegin;
4498   /* currently sctx is memory bleed */
4499   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
4500   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4501   /*
4502      This should work, but it doesn't
4503   sctx->ctx = ctx;
4504   mexMakeArrayPersistent(sctx->ctx);
4505   */
4506   sctx->ctx = mxDuplicateArray(ctx);
4507 
4508   ierr = TSSetIFunction(ts,NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
4509   PetscFunctionReturn(0);
4510 }
4511 
4512 #undef __FUNCT__
4513 #define __FUNCT__ "TSComputeJacobian_Matlab"
4514 /*
4515    TSComputeJacobian_Matlab - Calls the function that has been set with
4516                          TSSetJacobianMatlab().
4517 
4518    Collective on TS
4519 
4520    Input Parameters:
4521 +  ts - the TS context
4522 .  u - input vector
4523 .  A, B - the matrices
4524 -  ctx - user context
4525 
4526    Output Parameter:
4527 .  flag - structure of the matrix
4528 
4529    Level: developer
4530 
4531 .keywords: TS, nonlinear, compute, function
4532 
4533 .seealso: TSSetFunction(), TSGetFunction()
4534 @*/
4535 PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec u,Vec udot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
4536 {
4537   PetscErrorCode  ierr;
4538   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
4539   int             nlhs  = 2,nrhs = 9;
4540   mxArray         *plhs[2],*prhs[9];
4541   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
4542 
4543   PetscFunctionBegin;
4544   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4545   PetscValidHeaderSpecific(u,VEC_CLASSID,3);
4546 
4547   /* call Matlab function in ctx with arguments u and y */
4548 
4549   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
4550   ierr = PetscMemcpy(&lx,&u,sizeof(u));CHKERRQ(ierr);
4551   ierr = PetscMemcpy(&lxdot,&udot,sizeof(u));CHKERRQ(ierr);
4552   ierr = PetscMemcpy(&lA,A,sizeof(u));CHKERRQ(ierr);
4553   ierr = PetscMemcpy(&lB,B,sizeof(u));CHKERRQ(ierr);
4554 
4555   prhs[0] =  mxCreateDoubleScalar((double)ls);
4556   prhs[1] =  mxCreateDoubleScalar((double)time);
4557   prhs[2] =  mxCreateDoubleScalar((double)lx);
4558   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
4559   prhs[4] =  mxCreateDoubleScalar((double)shift);
4560   prhs[5] =  mxCreateDoubleScalar((double)lA);
4561   prhs[6] =  mxCreateDoubleScalar((double)lB);
4562   prhs[7] =  mxCreateString(sctx->funcname);
4563   prhs[8] =  sctx->ctx;
4564   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
4565   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4566   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
4567   mxDestroyArray(prhs[0]);
4568   mxDestroyArray(prhs[1]);
4569   mxDestroyArray(prhs[2]);
4570   mxDestroyArray(prhs[3]);
4571   mxDestroyArray(prhs[4]);
4572   mxDestroyArray(prhs[5]);
4573   mxDestroyArray(prhs[6]);
4574   mxDestroyArray(prhs[7]);
4575   mxDestroyArray(plhs[0]);
4576   mxDestroyArray(plhs[1]);
4577   PetscFunctionReturn(0);
4578 }
4579 
4580 
4581 #undef __FUNCT__
4582 #define __FUNCT__ "TSSetJacobianMatlab"
4583 /*
4584    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
4585    vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function
4586 
4587    Logically Collective on TS
4588 
4589    Input Parameters:
4590 +  ts - the TS context
4591 .  A,B - Jacobian matrices
4592 .  func - function evaluation routine
4593 -  ctx - user context
4594 
4595    Calling sequence of func:
4596 $    flag = func (TS ts,PetscReal time,Vec u,Vec udot,Mat A,Mat B,void *ctx);
4597 
4598 
4599    Level: developer
4600 
4601 .keywords: TS, nonlinear, set, function
4602 
4603 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
4604 */
4605 PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
4606 {
4607   PetscErrorCode  ierr;
4608   TSMatlabContext *sctx;
4609 
4610   PetscFunctionBegin;
4611   /* currently sctx is memory bleed */
4612   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
4613   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4614   /*
4615      This should work, but it doesn't
4616   sctx->ctx = ctx;
4617   mexMakeArrayPersistent(sctx->ctx);
4618   */
4619   sctx->ctx = mxDuplicateArray(ctx);
4620 
4621   ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
4622   PetscFunctionReturn(0);
4623 }
4624 
4625 #undef __FUNCT__
4626 #define __FUNCT__ "TSMonitor_Matlab"
4627 /*
4628    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
4629 
4630    Collective on TS
4631 
4632 .seealso: TSSetFunction(), TSGetFunction()
4633 @*/
4634 PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec u, void *ctx)
4635 {
4636   PetscErrorCode  ierr;
4637   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
4638   int             nlhs  = 1,nrhs = 6;
4639   mxArray         *plhs[1],*prhs[6];
4640   long long int   lx = 0,ls = 0;
4641 
4642   PetscFunctionBegin;
4643   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4644   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
4645 
4646   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
4647   ierr = PetscMemcpy(&lx,&u,sizeof(u));CHKERRQ(ierr);
4648 
4649   prhs[0] =  mxCreateDoubleScalar((double)ls);
4650   prhs[1] =  mxCreateDoubleScalar((double)it);
4651   prhs[2] =  mxCreateDoubleScalar((double)time);
4652   prhs[3] =  mxCreateDoubleScalar((double)lx);
4653   prhs[4] =  mxCreateString(sctx->funcname);
4654   prhs[5] =  sctx->ctx;
4655   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
4656   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4657   mxDestroyArray(prhs[0]);
4658   mxDestroyArray(prhs[1]);
4659   mxDestroyArray(prhs[2]);
4660   mxDestroyArray(prhs[3]);
4661   mxDestroyArray(prhs[4]);
4662   mxDestroyArray(plhs[0]);
4663   PetscFunctionReturn(0);
4664 }
4665 
4666 
4667 #undef __FUNCT__
4668 #define __FUNCT__ "TSMonitorSetMatlab"
4669 /*
4670    TSMonitorSetMatlab - Sets the monitor function from Matlab
4671 
4672    Level: developer
4673 
4674 .keywords: TS, nonlinear, set, function
4675 
4676 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
4677 */
4678 PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
4679 {
4680   PetscErrorCode  ierr;
4681   TSMatlabContext *sctx;
4682 
4683   PetscFunctionBegin;
4684   /* currently sctx is memory bleed */
4685   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
4686   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4687   /*
4688      This should work, but it doesn't
4689   sctx->ctx = ctx;
4690   mexMakeArrayPersistent(sctx->ctx);
4691   */
4692   sctx->ctx = mxDuplicateArray(ctx);
4693 
4694   ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,NULL);CHKERRQ(ierr);
4695   PetscFunctionReturn(0);
4696 }
4697 #endif
4698 
4699 
4700 
4701 #undef __FUNCT__
4702 #define __FUNCT__ "TSMonitorLGSolution"
4703 /*@C
4704    TSMonitorLGSolution - Monitors progress of the TS solvers by plotting each component of the solution vector
4705        in a time based line graph
4706 
4707    Collective on TS
4708 
4709    Input Parameters:
4710 +  ts - the TS context
4711 .  step - current time-step
4712 .  ptime - current time
4713 -  lg - a line graph object
4714 
4715    Level: intermediate
4716 
4717     Notes: each process in a parallel run displays its component solutions in a separate window
4718 
4719 .keywords: TS,  vector, monitor, view
4720 
4721 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4722 @*/
4723 PetscErrorCode  TSMonitorLGSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4724 {
4725   PetscErrorCode    ierr;
4726   TSMonitorLGCtx    ctx = (TSMonitorLGCtx)dummy;
4727   const PetscScalar *yy;
4728   PetscInt          dim;
4729 
4730   PetscFunctionBegin;
4731   if (!step) {
4732     PetscDrawAxis axis;
4733     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
4734     ierr = PetscDrawAxisSetLabels(axis,"Solution as function of time","Time","Solution");CHKERRQ(ierr);
4735     ierr = VecGetLocalSize(u,&dim);CHKERRQ(ierr);
4736     ierr = PetscDrawLGSetDimension(ctx->lg,dim);CHKERRQ(ierr);
4737     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
4738   }
4739   ierr = VecGetArrayRead(u,&yy);CHKERRQ(ierr);
4740 #if defined(PETSC_USE_COMPLEX)
4741   {
4742     PetscReal *yreal;
4743     PetscInt  i,n;
4744     ierr = VecGetLocalSize(u,&n);CHKERRQ(ierr);
4745     ierr = PetscMalloc(n*sizeof(PetscReal),&yreal);CHKERRQ(ierr);
4746     for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
4747     ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);CHKERRQ(ierr);
4748     ierr = PetscFree(yreal);CHKERRQ(ierr);
4749   }
4750 #else
4751   ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);CHKERRQ(ierr);
4752 #endif
4753   ierr = VecRestoreArrayRead(u,&yy);CHKERRQ(ierr);
4754   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
4755     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
4756   }
4757   PetscFunctionReturn(0);
4758 }
4759 
4760 #undef __FUNCT__
4761 #define __FUNCT__ "TSMonitorLGError"
4762 /*@C
4763    TSMonitorLGError - Monitors progress of the TS solvers by plotting each component of the solution vector
4764        in a time based line graph
4765 
4766    Collective on TS
4767 
4768    Input Parameters:
4769 +  ts - the TS context
4770 .  step - current time-step
4771 .  ptime - current time
4772 -  lg - a line graph object
4773 
4774    Level: intermediate
4775 
4776    Notes:
4777    Only for sequential solves.
4778 
4779    The user must provide the solution using TSSetSolutionFunction() to use this monitor.
4780 
4781    Options Database Keys:
4782 .  -ts_monitor_lg_error - create a graphical monitor of error history
4783 
4784 .keywords: TS,  vector, monitor, view
4785 
4786 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
4787 @*/
4788 PetscErrorCode  TSMonitorLGError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4789 {
4790   PetscErrorCode    ierr;
4791   TSMonitorLGCtx    ctx = (TSMonitorLGCtx)dummy;
4792   const PetscScalar *yy;
4793   Vec               y;
4794   PetscInt          dim;
4795 
4796   PetscFunctionBegin;
4797   if (!step) {
4798     PetscDrawAxis axis;
4799     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
4800     ierr = PetscDrawAxisSetLabels(axis,"Error in solution as function of time","Time","Solution");CHKERRQ(ierr);
4801     ierr = VecGetLocalSize(u,&dim);CHKERRQ(ierr);
4802     ierr = PetscDrawLGSetDimension(ctx->lg,dim);CHKERRQ(ierr);
4803     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
4804   }
4805   ierr = VecDuplicate(u,&y);CHKERRQ(ierr);
4806   ierr = TSComputeSolutionFunction(ts,ptime,y);CHKERRQ(ierr);
4807   ierr = VecAXPY(y,-1.0,u);CHKERRQ(ierr);
4808   ierr = VecGetArrayRead(y,&yy);CHKERRQ(ierr);
4809 #if defined(PETSC_USE_COMPLEX)
4810   {
4811     PetscReal *yreal;
4812     PetscInt  i,n;
4813     ierr = VecGetLocalSize(y,&n);CHKERRQ(ierr);
4814     ierr = PetscMalloc(n*sizeof(PetscReal),&yreal);CHKERRQ(ierr);
4815     for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
4816     ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);CHKERRQ(ierr);
4817     ierr = PetscFree(yreal);CHKERRQ(ierr);
4818   }
4819 #else
4820   ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);CHKERRQ(ierr);
4821 #endif
4822   ierr = VecRestoreArrayRead(y,&yy);CHKERRQ(ierr);
4823   ierr = VecDestroy(&y);CHKERRQ(ierr);
4824   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
4825     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
4826   }
4827   PetscFunctionReturn(0);
4828 }
4829 
4830 #undef __FUNCT__
4831 #define __FUNCT__ "TSMonitorLGSNESIterations"
4832 PetscErrorCode TSMonitorLGSNESIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
4833 {
4834   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
4835   PetscReal      x   = ptime,y;
4836   PetscErrorCode ierr;
4837   PetscInt       its;
4838 
4839   PetscFunctionBegin;
4840   if (!n) {
4841     PetscDrawAxis axis;
4842 
4843     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
4844     ierr = PetscDrawAxisSetLabels(axis,"Nonlinear iterations as function of time","Time","SNES Iterations");CHKERRQ(ierr);
4845     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
4846 
4847     ctx->snes_its = 0;
4848   }
4849   ierr = TSGetSNESIterations(ts,&its);CHKERRQ(ierr);
4850   y    = its - ctx->snes_its;
4851   ierr = PetscDrawLGAddPoint(ctx->lg,&x,&y);CHKERRQ(ierr);
4852   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
4853     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
4854   }
4855   ctx->snes_its = its;
4856   PetscFunctionReturn(0);
4857 }
4858 
4859 #undef __FUNCT__
4860 #define __FUNCT__ "TSMonitorLGKSPIterations"
4861 PetscErrorCode TSMonitorLGKSPIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
4862 {
4863   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
4864   PetscReal      x   = ptime,y;
4865   PetscErrorCode ierr;
4866   PetscInt       its;
4867 
4868   PetscFunctionBegin;
4869   if (!n) {
4870     PetscDrawAxis axis;
4871 
4872     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
4873     ierr = PetscDrawAxisSetLabels(axis,"Linear iterations as function of time","Time","KSP Iterations");CHKERRQ(ierr);
4874     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
4875 
4876     ctx->ksp_its = 0;
4877   }
4878   ierr = TSGetKSPIterations(ts,&its);CHKERRQ(ierr);
4879   y    = its - ctx->ksp_its;
4880   ierr = PetscDrawLGAddPoint(ctx->lg,&x,&y);CHKERRQ(ierr);
4881   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
4882     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
4883   }
4884   ctx->ksp_its = its;
4885   PetscFunctionReturn(0);
4886 }
4887 
4888 #undef __FUNCT__
4889 #define __FUNCT__ "TSComputeLinearStability"
4890 /*@
4891    TSComputeLinearStability - computes the linear stability function at a point
4892 
4893    Collective on TS and Vec
4894 
4895    Input Parameters:
4896 +  ts - the TS context
4897 -  xr,xi - real and imaginary part of input arguments
4898 
4899    Output Parameters:
4900 .  yr,yi - real and imaginary part of function value
4901 
4902    Level: developer
4903 
4904 .keywords: TS, compute
4905 
4906 .seealso: TSSetRHSFunction(), TSComputeIFunction()
4907 @*/
4908 PetscErrorCode TSComputeLinearStability(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
4909 {
4910   PetscErrorCode ierr;
4911 
4912   PetscFunctionBegin;
4913   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4914   if (!ts->ops->linearstability) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Linearized stability function not provided for this method");
4915   ierr = (*ts->ops->linearstability)(ts,xr,xi,yr,yi);CHKERRQ(ierr);
4916   PetscFunctionReturn(0);
4917 }
4918