xref: /petsc/src/ts/interface/ts.c (revision c094ef4021e955ef5f85f7d8a1bbc6ed64ba7621)
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 struct _n_TSMonitorDrawCtx {
15   PetscViewer   viewer;
16   PetscDrawAxis axis;
17   Vec           initialsolution;
18   PetscBool     showinitial;
19   PetscInt      howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
20   PetscBool     showtimestepandtime;
21   int           color;
22 };
23 
24 #undef __FUNCT__
25 #define __FUNCT__ "TSMonitorSetFromOptions"
26 /*@C
27    TSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
28 
29    Collective on TS
30 
31    Input Parameters:
32 +  ts - TS object you wish to monitor
33 .  name - the monitor type one is seeking
34 .  help - message indicating what monitoring is done
35 .  manual - manual page for the monitor
36 .  monitor - the monitor function
37 -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects
38 
39    Level: developer
40 
41 .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
42           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
43           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
44           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
45           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
46           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
47           PetscOptionsFList(), PetscOptionsEList()
48 @*/
49 PetscErrorCode  TSMonitorSetFromOptions(TS ts,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),PetscErrorCode (*monitorsetup)(TS,PetscViewer))
50 {
51   PetscErrorCode    ierr;
52   PetscViewer       viewer;
53   PetscViewerFormat format;
54   PetscBool         flg;
55 
56   PetscFunctionBegin;
57   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
58   if (flg) {
59     ierr = PetscViewerSetFormat(viewer,format);CHKERRQ(ierr);
60     if (monitorsetup) {
61       ierr = (*monitorsetup)(ts,viewer);CHKERRQ(ierr);
62     }
63     ierr = TSMonitorSet(ts,monitor,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
64   }
65   PetscFunctionReturn(0);
66 }
67 
68 #undef __FUNCT__
69 #define __FUNCT__ "TSAdjointMonitorSetFromOptions"
70 /*@C
71    TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
72 
73    Collective on TS
74 
75    Input Parameters:
76 +  ts - TS object you wish to monitor
77 .  name - the monitor type one is seeking
78 .  help - message indicating what monitoring is done
79 .  manual - manual page for the monitor
80 .  monitor - the monitor function
81 -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects
82 
83    Level: developer
84 
85 .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
86           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
87           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
88           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
89           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
90           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
91           PetscOptionsFList(), PetscOptionsEList()
92 @*/
93 PetscErrorCode  TSAdjointMonitorSetFromOptions(TS ts,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),PetscErrorCode (*monitorsetup)(TS,PetscViewer))
94 {
95   PetscErrorCode    ierr;
96   PetscViewer       viewer;
97   PetscViewerFormat format;
98   PetscBool         flg;
99 
100   PetscFunctionBegin;
101   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
102   if (flg) {
103     ierr = PetscViewerSetFormat(viewer,format);CHKERRQ(ierr);
104     if (monitorsetup) {
105       ierr = (*monitorsetup)(ts,viewer);CHKERRQ(ierr);
106     }
107     ierr = TSAdjointMonitorSet(ts,monitor,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
108   }
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "TSSetFromOptions"
114 /*@
115    TSSetFromOptions - Sets various TS parameters from user options.
116 
117    Collective on TS
118 
119    Input Parameter:
120 .  ts - the TS context obtained from TSCreate()
121 
122    Options Database Keys:
123 +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP
124 .  -ts_save_trajectory - checkpoint the solution at each time-step
125 .  -ts_max_steps <maxsteps> - maximum number of time-steps to take
126 .  -ts_final_time <time> - maximum time to compute to
127 .  -ts_dt <dt> - initial time step
128 .  -ts_exact_final_time <stepover,interpolate,matchstep> whether to stop at the exact given final time and how to compute the solution at that ti,e
129 .  -ts_max_snes_failures <maxfailures> - Maximum number of nonlinear solve failures allowed
130 .  -ts_max_reject <maxrejects> - Maximum number of step rejections before step fails
131 .  -ts_error_if_step_fails <true,false> - Error if no step succeeds
132 .  -ts_rtol <rtol> - relative tolerance for local truncation error
133 .  -ts_atol <atol> Absolute tolerance for local truncation error
134 .  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
135 .  -ts_fd_color - Use finite differences with coloring to compute IJacobian
136 .  -ts_monitor - print information at each timestep
137 .  -ts_monitor_lg_timestep - Monitor timestep size graphically
138 .  -ts_monitor_lg_solution - Monitor solution graphically
139 .  -ts_monitor_lg_error - Monitor error graphically
140 .  -ts_monitor_lg_snes_iterations - Monitor number nonlinear iterations for each timestep graphically
141 .  -ts_monitor_lg_ksp_iterations - Monitor number nonlinear iterations for each timestep graphically
142 .  -ts_monitor_sp_eig - Monitor eigenvalues of linearized operator graphically
143 .  -ts_monitor_draw_solution - Monitor solution graphically
144 .  -ts_monitor_draw_solution_phase  <xleft,yleft,xright,yright> - Monitor solution graphically with phase diagram, requires problem with exactly 2 degrees of freedom
145 .  -ts_monitor_draw_error - Monitor error graphically, requires use to have provided TSSetSolutionFunction()
146 .  -ts_monitor_solution [ascii binary draw][:filename][:viewerformat] - monitors the solution at each timestep
147 .  -ts_monitor_solution_vtk <filename.vts> - Save each time step to a binary file, use filename-%%03D.vts
148 .  -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
149 .  -ts_adjoint_monitor - print information at each adjoint time step
150 -  -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
151 
152    Developer Note: We should unify all the -ts_monitor options in the way that -xxx_view has been unified
153 
154    Level: beginner
155 
156 .keywords: TS, timestep, set, options, database
157 
158 .seealso: TSGetType()
159 @*/
160 PetscErrorCode  TSSetFromOptions(TS ts)
161 {
162   PetscBool              opt,flg,tflg;
163   PetscErrorCode         ierr;
164   char                   monfilename[PETSC_MAX_PATH_LEN];
165   SNES                   snes;
166   TSAdapt                adapt;
167   PetscReal              time_step;
168   TSExactFinalTimeOption eftopt;
169   char                   dir[16];
170   const char             *defaultType;
171   char                   typeName[256];
172 
173   PetscFunctionBegin;
174   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
175   ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr);
176   if (((PetscObject)ts)->type_name) defaultType = ((PetscObject)ts)->type_name;
177   else defaultType = TSEULER;
178 
179   ierr = TSRegisterAll();CHKERRQ(ierr);
180   ierr = PetscOptionsFList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
181   if (opt) {
182     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
183   } else {
184     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
185   }
186 
187   /* Handle generic TS options */
188   if (ts->trajectory) tflg = PETSC_TRUE;
189   else tflg = PETSC_FALSE;
190   ierr = PetscOptionsBool("-ts_save_trajectory","Save the solution at each timestep","TSSetSaveTrajectory",tflg,&tflg,NULL);CHKERRQ(ierr);
191   if (tflg) {ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);}
192   if (ts->adjoint_solve) tflg = PETSC_TRUE;
193   else tflg = PETSC_FALSE;
194   ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&flg);CHKERRQ(ierr);
195   if (flg) {
196     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
197     ts->adjoint_solve = tflg;
198   }
199   ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,NULL);CHKERRQ(ierr);
200   ierr = PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,NULL);CHKERRQ(ierr);
201   ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,NULL);CHKERRQ(ierr);
202   ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&time_step,&flg);CHKERRQ(ierr);
203   if (flg) {
204     ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
205   }
206   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);
207   if (flg) {ierr = TSSetExactFinalTime(ts,eftopt);CHKERRQ(ierr);}
208   ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,NULL);CHKERRQ(ierr);
209   ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,NULL);CHKERRQ(ierr);
210   ierr = PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,NULL);CHKERRQ(ierr);
211   ierr = PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,NULL);CHKERRQ(ierr);
212   ierr = PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,NULL);CHKERRQ(ierr);
213 
214 #if defined(PETSC_HAVE_SAWS)
215   {
216   PetscBool set;
217   flg  = PETSC_FALSE;
218   ierr = PetscOptionsBool("-ts_saws_block","Block for SAWs memory snooper at end of TSSolve","PetscObjectSAWsBlock",((PetscObject)ts)->amspublishblock,&flg,&set);CHKERRQ(ierr);
219   if (set) {
220     ierr = PetscObjectSAWsSetBlock((PetscObject)ts,flg);CHKERRQ(ierr);
221   }
222   }
223 #endif
224 
225   /* Monitor options */
226   ierr = TSMonitorSetFromOptions(ts,"-ts_monitor","Monitor timestep size","TSMonitorDefault",TSMonitorDefault,NULL);CHKERRQ(ierr);
227   ierr = TSMonitorSetFromOptions(ts,"-ts_monitor_solution","View the solution at each timestep","TSMonitorSolution",TSMonitorSolution,NULL);CHKERRQ(ierr);
228   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr);
229 
230   ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
231   if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);}
232 
233   ierr = PetscOptionsName("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",&opt);CHKERRQ(ierr);
234   if (opt) {
235     TSMonitorLGCtx ctx;
236     PetscInt       howoften = 1;
237 
238     ierr = PetscOptionsInt("-ts_monitor_lg_timestep","Monitor timestep size graphically","TSMonitorLGTimeStep",howoften,&howoften,NULL);CHKERRQ(ierr);
239     ierr = TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
240     ierr = TSMonitorSet(ts,TSMonitorLGTimeStep,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
241   }
242   ierr = PetscOptionsName("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",&opt);CHKERRQ(ierr);
243   if (opt) {
244     TSMonitorLGCtx ctx;
245     PetscInt       howoften = 1;
246 
247     ierr = PetscOptionsInt("-ts_monitor_lg_solution","Monitor solution graphically","TSMonitorLGSolution",howoften,&howoften,NULL);CHKERRQ(ierr);
248     ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
249     ierr = TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
250   }
251   ierr = PetscOptionsName("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",&opt);CHKERRQ(ierr);
252   if (opt) {
253     TSMonitorLGCtx ctx;
254     PetscInt       howoften = 1;
255 
256     ierr = PetscOptionsInt("-ts_monitor_lg_error","Monitor error graphically","TSMonitorLGError",howoften,&howoften,NULL);CHKERRQ(ierr);
257     ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
258     ierr = TSMonitorSet(ts,TSMonitorLGError,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
259   }
260   ierr = PetscOptionsName("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",&opt);CHKERRQ(ierr);
261   if (opt) {
262     TSMonitorLGCtx ctx;
263     PetscInt       howoften = 1;
264 
265     ierr = PetscOptionsInt("-ts_monitor_lg_snes_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGSNESIterations",howoften,&howoften,NULL);CHKERRQ(ierr);
266     ierr = TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
267     ierr = TSMonitorSet(ts,TSMonitorLGSNESIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
268   }
269   ierr = PetscOptionsName("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",&opt);CHKERRQ(ierr);
270   if (opt) {
271     TSMonitorLGCtx ctx;
272     PetscInt       howoften = 1;
273 
274     ierr = PetscOptionsInt("-ts_monitor_lg_ksp_iterations","Monitor number nonlinear iterations for each timestep graphically","TSMonitorLGKSPIterations",howoften,&howoften,NULL);CHKERRQ(ierr);
275     ierr = TSMonitorLGCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
276     ierr = TSMonitorSet(ts,TSMonitorLGKSPIterations,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr);
277   }
278   ierr = PetscOptionsName("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",&opt);CHKERRQ(ierr);
279   if (opt) {
280     TSMonitorSPEigCtx ctx;
281     PetscInt          howoften = 1;
282 
283     ierr = PetscOptionsInt("-ts_monitor_sp_eig","Monitor eigenvalues of linearized operator graphically","TSMonitorSPEig",howoften,&howoften,NULL);CHKERRQ(ierr);
284     ierr = TSMonitorSPEigCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
285     ierr = TSMonitorSet(ts,TSMonitorSPEig,ctx,(PetscErrorCode (*)(void**))TSMonitorSPEigCtxDestroy);CHKERRQ(ierr);
286   }
287   opt  = PETSC_FALSE;
288   ierr = PetscOptionsName("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",&opt);CHKERRQ(ierr);
289   if (opt) {
290     TSMonitorDrawCtx ctx;
291     PetscInt         howoften = 1;
292 
293     ierr = PetscOptionsInt("-ts_monitor_draw_solution","Monitor solution graphically","TSMonitorDrawSolution",howoften,&howoften,NULL);CHKERRQ(ierr);
294     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
295     ierr = TSMonitorSet(ts,TSMonitorDrawSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
296   }
297   opt  = PETSC_FALSE;
298   ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr);
299   if (opt) {
300     TSMonitorDrawCtx ctx;
301     PetscInt         howoften = 1;
302 
303     ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr);
304     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
305     ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
306   }
307   opt  = PETSC_FALSE;
308   ierr = PetscOptionsName("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",&opt);CHKERRQ(ierr);
309   if (opt) {
310     TSMonitorDrawCtx ctx;
311     PetscReal        bounds[4];
312     PetscInt         n = 4;
313     PetscDraw        draw;
314 
315     ierr = PetscOptionsRealArray("-ts_monitor_draw_solution_phase","Monitor solution graphically","TSMonitorDrawSolutionPhase",bounds,&n,NULL);CHKERRQ(ierr);
316     if (n != 4) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Must provide bounding box of phase field");
317     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx);CHKERRQ(ierr);
318     ierr = PetscViewerDrawGetDraw(ctx->viewer,0,&draw);CHKERRQ(ierr);
319     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
320     ierr = PetscDrawAxisCreate(draw,&ctx->axis);CHKERRQ(ierr);
321     ierr = PetscDrawAxisSetLimits(ctx->axis,bounds[0],bounds[2],bounds[1],bounds[3]);CHKERRQ(ierr);
322     ierr = PetscDrawAxisSetLabels(ctx->axis,"Phase Diagram","Variable 1","Variable 2");CHKERRQ(ierr);
323     ierr = PetscDrawAxisDraw(ctx->axis);CHKERRQ(ierr);
324     /* ierr = PetscDrawSetCoordinates(draw,bounds[0],bounds[1],bounds[2],bounds[3]);CHKERRQ(ierr); */
325     ierr = TSMonitorSet(ts,TSMonitorDrawSolutionPhase,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
326   }
327   opt  = PETSC_FALSE;
328   ierr = PetscOptionsName("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",&opt);CHKERRQ(ierr);
329   if (opt) {
330     TSMonitorDrawCtx ctx;
331     PetscInt         howoften = 1;
332 
333     ierr = PetscOptionsInt("-ts_monitor_draw_error","Monitor error graphically","TSMonitorDrawError",howoften,&howoften,NULL);CHKERRQ(ierr);
334     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&ctx);CHKERRQ(ierr);
335     ierr = TSMonitorSet(ts,TSMonitorDrawError,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
336   }
337 
338   opt  = PETSC_FALSE;
339   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);
340   if (flg) {
341     const char *ptr,*ptr2;
342     char       *filetemplate;
343     if (!monfilename[0]) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
344     /* Do some cursory validation of the input. */
345     ierr = PetscStrstr(monfilename,"%",(char**)&ptr);CHKERRQ(ierr);
346     if (!ptr) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
347     for (ptr++; ptr && *ptr; ptr++) {
348       ierr = PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);CHKERRQ(ierr);
349       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");
350       if (ptr2) break;
351     }
352     ierr = PetscStrallocpy(monfilename,&filetemplate);CHKERRQ(ierr);
353     ierr = TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);CHKERRQ(ierr);
354   }
355 
356   ierr = PetscOptionsString("-ts_monitor_dmda_ray","Display a ray of the solution","None","y=0",dir,16,&flg);CHKERRQ(ierr);
357   if (flg) {
358     TSMonitorDMDARayCtx *rayctx;
359     int                  ray = 0;
360     DMDADirection        ddir;
361     DM                   da;
362     PetscMPIInt          rank;
363 
364     if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
365     if (dir[0] == 'x') ddir = DMDA_X;
366     else if (dir[0] == 'y') ddir = DMDA_Y;
367     else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Unknown ray %s",dir);
368     sscanf(dir+2,"%d",&ray);
369 
370     ierr = PetscInfo2(((PetscObject)ts),"Displaying DMDA ray %c = %D\n",dir[0],ray);CHKERRQ(ierr);
371     ierr = PetscNew(&rayctx);CHKERRQ(ierr);
372     ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
373     ierr = DMDAGetRay(da,ddir,ray,&rayctx->ray,&rayctx->scatter);CHKERRQ(ierr);
374     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
375     if (!rank) {
376       ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,600,300,&rayctx->viewer);CHKERRQ(ierr);
377     }
378     rayctx->lgctx = NULL;
379     ierr = TSMonitorSet(ts,TSMonitorDMDARay,rayctx,TSMonitorDMDARayDestroy);CHKERRQ(ierr);
380   }
381   ierr = PetscOptionsString("-ts_monitor_lg_dmda_ray","Display a ray of the solution","None","x=0",dir,16,&flg);CHKERRQ(ierr);
382   if (flg) {
383     TSMonitorDMDARayCtx *rayctx;
384     int                 ray = 0;
385     DMDADirection       ddir;
386     DM                  da;
387     PetscInt            howoften = 1;
388 
389     if (dir[1] != '=') SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Malformed ray %s", dir);
390     if      (dir[0] == 'x') ddir = DMDA_X;
391     else if (dir[0] == 'y') ddir = DMDA_Y;
392     else SETERRQ1(PetscObjectComm((PetscObject) ts), PETSC_ERR_ARG_WRONG, "Unknown ray direction %s", dir);
393     sscanf(dir+2, "%d", &ray);
394 
395     ierr = PetscInfo2(((PetscObject) ts),"Displaying LG DMDA ray %c = %D\n", dir[0], ray);CHKERRQ(ierr);
396     ierr = PetscNew(&rayctx);CHKERRQ(ierr);
397     ierr = TSGetDM(ts, &da);CHKERRQ(ierr);
398     ierr = DMDAGetRay(da, ddir, ray, &rayctx->ray, &rayctx->scatter);CHKERRQ(ierr);
399     ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,howoften,&rayctx->lgctx);CHKERRQ(ierr);
400     ierr = TSMonitorSet(ts, TSMonitorLGDMDARay, rayctx, TSMonitorDMDARayDestroy);CHKERRQ(ierr);
401   }
402 
403   ierr = PetscOptionsName("-ts_monitor_envelope","Monitor maximum and minimum value of each component of the solution","TSMonitorEnvelope",&opt);CHKERRQ(ierr);
404   if (opt) {
405     TSMonitorEnvelopeCtx ctx;
406 
407     ierr = TSMonitorEnvelopeCtxCreate(ts,&ctx);CHKERRQ(ierr);
408     ierr = TSMonitorSet(ts,TSMonitorEnvelope,ctx,(PetscErrorCode (*)(void**))TSMonitorEnvelopeCtxDestroy);CHKERRQ(ierr);
409   }
410 
411   flg  = PETSC_FALSE;
412   ierr = PetscOptionsBool("-ts_fd_color", "Use finite differences with coloring to compute IJacobian", "TSComputeJacobianDefaultColor", flg, &flg, NULL);CHKERRQ(ierr);
413   if (flg) {
414     DM   dm;
415     DMTS tdm;
416 
417     ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
418     ierr = DMGetDMTS(dm, &tdm);CHKERRQ(ierr);
419     tdm->ijacobianctx = NULL;
420     ierr = TSSetIJacobian(ts, NULL, NULL, TSComputeIJacobianDefaultColor, 0);CHKERRQ(ierr);
421     ierr = PetscInfo(ts, "Setting default finite difference coloring Jacobian matrix\n");CHKERRQ(ierr);
422   }
423 
424   /*
425      This code is all wrong. One is creating objects inside the TSSetFromOptions() so if run with the options gui
426      will bleed memory. Also one is using a PetscOptionsBegin() inside a PetscOptionsBegin()
427   */
428   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
429   ierr = TSAdaptSetFromOptions(PetscOptionsObject,adapt);CHKERRQ(ierr);
430 
431     /* Handle specific TS options */
432   if (ts->ops->setfromoptions) {
433     ierr = (*ts->ops->setfromoptions)(PetscOptionsObject,ts);CHKERRQ(ierr);
434   }
435 
436   /* process any options handlers added with PetscObjectAddOptionsHandler() */
437   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)ts);CHKERRQ(ierr);
438   ierr = PetscOptionsEnd();CHKERRQ(ierr);
439 
440   if (ts->trajectory) {
441     ierr = TSTrajectorySetFromOptions(ts->trajectory);CHKERRQ(ierr);
442   }
443 
444   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
445   if (snes) {
446     if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
447     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
448   }
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "TSSetSaveTrajectory"
454 /*@
455    TSSetSaveTrajectory - Causes the TS to save its solutions as it iterates forward in time in a TSTrajectory object
456 
457    Collective on TS
458 
459    Input Parameters:
460 .  ts - the TS context obtained from TSCreate()
461 
462 
463    Level: intermediate
464 
465 .seealso: TSGetTrajectory(), TSAdjointSolve()
466 
467 .keywords: TS, set, checkpoint,
468 @*/
469 PetscErrorCode  TSSetSaveTrajectory(TS ts)
470 {
471   PetscErrorCode ierr;
472 
473   PetscFunctionBegin;
474   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
475   if (!ts->trajectory) {
476     ierr = TSTrajectoryCreate(PetscObjectComm((PetscObject)ts),&ts->trajectory);CHKERRQ(ierr);
477     ierr = TSTrajectorySetType(ts->trajectory,TSTRAJECTORYBASIC);CHKERRQ(ierr);
478     ierr = TSTrajectorySetFromOptions(ts->trajectory);CHKERRQ(ierr);
479   }
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "TSComputeRHSJacobian"
485 /*@
486    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
487       set with TSSetRHSJacobian().
488 
489    Collective on TS and Vec
490 
491    Input Parameters:
492 +  ts - the TS context
493 .  t - current timestep
494 -  U - input vector
495 
496    Output Parameters:
497 +  A - Jacobian matrix
498 .  B - optional preconditioning matrix
499 -  flag - flag indicating matrix structure
500 
501    Notes:
502    Most users should not need to explicitly call this routine, as it
503    is used internally within the nonlinear solvers.
504 
505    See KSPSetOperators() for important information about setting the
506    flag parameter.
507 
508    Level: developer
509 
510 .keywords: SNES, compute, Jacobian, matrix
511 
512 .seealso:  TSSetRHSJacobian(), KSPSetOperators()
513 @*/
514 PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B)
515 {
516   PetscErrorCode ierr;
517   PetscObjectState Ustate;
518   DM             dm;
519   DMTS           tsdm;
520   TSRHSJacobian  rhsjacobianfunc;
521   void           *ctx;
522   TSIJacobian    ijacobianfunc;
523   TSRHSFunction  rhsfunction;
524 
525   PetscFunctionBegin;
526   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
527   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
528   PetscCheckSameComm(ts,1,U,3);
529   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
530   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
531   ierr = DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);CHKERRQ(ierr);
532   ierr = DMTSGetIJacobian(dm,&ijacobianfunc,NULL);CHKERRQ(ierr);
533   ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctx);CHKERRQ(ierr);
534   ierr = PetscObjectStateGet((PetscObject)U,&Ustate);CHKERRQ(ierr);
535   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == U && ts->rhsjacobian.Xstate == Ustate)) && (rhsfunction != TSComputeRHSFunctionLinear)) {
536     PetscFunctionReturn(0);
537   }
538 
539   if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
540 
541   if (ts->rhsjacobian.reuse) {
542     ierr = MatShift(A,-ts->rhsjacobian.shift);CHKERRQ(ierr);
543     ierr = MatScale(A,1./ts->rhsjacobian.scale);CHKERRQ(ierr);
544     if (A != B) {
545       ierr = MatShift(B,-ts->rhsjacobian.shift);CHKERRQ(ierr);
546       ierr = MatScale(B,1./ts->rhsjacobian.scale);CHKERRQ(ierr);
547     }
548     ts->rhsjacobian.shift = 0;
549     ts->rhsjacobian.scale = 1.;
550   }
551 
552   if (rhsjacobianfunc) {
553     PetscBool missing;
554     ierr = PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);CHKERRQ(ierr);
555     PetscStackPush("TS user Jacobian function");
556     ierr = (*rhsjacobianfunc)(ts,t,U,A,B,ctx);CHKERRQ(ierr);
557     PetscStackPop;
558     ierr = PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);CHKERRQ(ierr);
559     if (A) {
560       ierr = MatMissingDiagonal(A,&missing,NULL);CHKERRQ(ierr);
561       if (missing) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Amat passed to TSSetRHSJacobian() must have all diagonal entries set, if they are zero you must still set them with a zero value");
562     }
563     if (B && B != A) {
564       ierr = MatMissingDiagonal(B,&missing,NULL);CHKERRQ(ierr);
565       if (missing) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Bmat passed to TSSetRHSJacobian() must have all diagonal entries set, if they are zero you must still set them with a zero value");
566     }
567   } else {
568     ierr = MatZeroEntries(A);CHKERRQ(ierr);
569     if (A != B) {ierr = MatZeroEntries(B);CHKERRQ(ierr);}
570   }
571   ts->rhsjacobian.time       = t;
572   ts->rhsjacobian.X          = U;
573   ierr                       = PetscObjectStateGet((PetscObject)U,&ts->rhsjacobian.Xstate);CHKERRQ(ierr);
574   PetscFunctionReturn(0);
575 }
576 
577 #undef __FUNCT__
578 #define __FUNCT__ "TSComputeRHSFunction"
579 /*@
580    TSComputeRHSFunction - Evaluates the right-hand-side function.
581 
582    Collective on TS and Vec
583 
584    Input Parameters:
585 +  ts - the TS context
586 .  t - current time
587 -  U - state vector
588 
589    Output Parameter:
590 .  y - right hand side
591 
592    Note:
593    Most users should not need to explicitly call this routine, as it
594    is used internally within the nonlinear solvers.
595 
596    Level: developer
597 
598 .keywords: TS, compute
599 
600 .seealso: TSSetRHSFunction(), TSComputeIFunction()
601 @*/
602 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec U,Vec y)
603 {
604   PetscErrorCode ierr;
605   TSRHSFunction  rhsfunction;
606   TSIFunction    ifunction;
607   void           *ctx;
608   DM             dm;
609 
610   PetscFunctionBegin;
611   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
612   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
613   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
614   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
615   ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctx);CHKERRQ(ierr);
616   ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr);
617 
618   if (!rhsfunction && !ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
619 
620   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,y,0);CHKERRQ(ierr);
621   if (rhsfunction) {
622     PetscStackPush("TS user right-hand-side function");
623     ierr = (*rhsfunction)(ts,t,U,y,ctx);CHKERRQ(ierr);
624     PetscStackPop;
625   } else {
626     ierr = VecZeroEntries(y);CHKERRQ(ierr);
627   }
628 
629   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,y,0);CHKERRQ(ierr);
630   PetscFunctionReturn(0);
631 }
632 
633 #undef __FUNCT__
634 #define __FUNCT__ "TSComputeSolutionFunction"
635 /*@
636    TSComputeSolutionFunction - Evaluates the solution function.
637 
638    Collective on TS and Vec
639 
640    Input Parameters:
641 +  ts - the TS context
642 -  t - current time
643 
644    Output Parameter:
645 .  U - the solution
646 
647    Note:
648    Most users should not need to explicitly call this routine, as it
649    is used internally within the nonlinear solvers.
650 
651    Level: developer
652 
653 .keywords: TS, compute
654 
655 .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
656 @*/
657 PetscErrorCode TSComputeSolutionFunction(TS ts,PetscReal t,Vec U)
658 {
659   PetscErrorCode     ierr;
660   TSSolutionFunction solutionfunction;
661   void               *ctx;
662   DM                 dm;
663 
664   PetscFunctionBegin;
665   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
666   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
667   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
668   ierr = DMTSGetSolutionFunction(dm,&solutionfunction,&ctx);CHKERRQ(ierr);
669 
670   if (solutionfunction) {
671     PetscStackPush("TS user solution function");
672     ierr = (*solutionfunction)(ts,t,U,ctx);CHKERRQ(ierr);
673     PetscStackPop;
674   }
675   PetscFunctionReturn(0);
676 }
677 #undef __FUNCT__
678 #define __FUNCT__ "TSComputeForcingFunction"
679 /*@
680    TSComputeForcingFunction - Evaluates the forcing function.
681 
682    Collective on TS and Vec
683 
684    Input Parameters:
685 +  ts - the TS context
686 -  t - current time
687 
688    Output Parameter:
689 .  U - the function value
690 
691    Note:
692    Most users should not need to explicitly call this routine, as it
693    is used internally within the nonlinear solvers.
694 
695    Level: developer
696 
697 .keywords: TS, compute
698 
699 .seealso: TSSetSolutionFunction(), TSSetRHSFunction(), TSComputeIFunction()
700 @*/
701 PetscErrorCode TSComputeForcingFunction(TS ts,PetscReal t,Vec U)
702 {
703   PetscErrorCode     ierr, (*forcing)(TS,PetscReal,Vec,void*);
704   void               *ctx;
705   DM                 dm;
706 
707   PetscFunctionBegin;
708   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
709   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
710   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
711   ierr = DMTSGetForcingFunction(dm,&forcing,&ctx);CHKERRQ(ierr);
712 
713   if (forcing) {
714     PetscStackPush("TS user forcing function");
715     ierr = (*forcing)(ts,t,U,ctx);CHKERRQ(ierr);
716     PetscStackPop;
717   }
718   PetscFunctionReturn(0);
719 }
720 
721 #undef __FUNCT__
722 #define __FUNCT__ "TSGetRHSVec_Private"
723 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
724 {
725   Vec            F;
726   PetscErrorCode ierr;
727 
728   PetscFunctionBegin;
729   *Frhs = NULL;
730   ierr  = TSGetIFunction(ts,&F,NULL,NULL);CHKERRQ(ierr);
731   if (!ts->Frhs) {
732     ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr);
733   }
734   *Frhs = ts->Frhs;
735   PetscFunctionReturn(0);
736 }
737 
738 #undef __FUNCT__
739 #define __FUNCT__ "TSGetRHSMats_Private"
740 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
741 {
742   Mat            A,B;
743   PetscErrorCode ierr;
744 
745   PetscFunctionBegin;
746   if (Arhs) *Arhs = NULL;
747   if (Brhs) *Brhs = NULL;
748   ierr = TSGetIJacobian(ts,&A,&B,NULL,NULL);CHKERRQ(ierr);
749   if (Arhs) {
750     if (!ts->Arhs) {
751       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
752     }
753     *Arhs = ts->Arhs;
754   }
755   if (Brhs) {
756     if (!ts->Brhs) {
757       if (A != B) {
758         ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
759       } else {
760         ts->Brhs = ts->Arhs;
761         ierr = PetscObjectReference((PetscObject)ts->Arhs);CHKERRQ(ierr);
762       }
763     }
764     *Brhs = ts->Brhs;
765   }
766   PetscFunctionReturn(0);
767 }
768 
769 #undef __FUNCT__
770 #define __FUNCT__ "TSComputeIFunction"
771 /*@
772    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,U,Udot)=0
773 
774    Collective on TS and Vec
775 
776    Input Parameters:
777 +  ts - the TS context
778 .  t - current time
779 .  U - state vector
780 .  Udot - time derivative of state vector
781 -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
782 
783    Output Parameter:
784 .  Y - right hand side
785 
786    Note:
787    Most users should not need to explicitly call this routine, as it
788    is used internally within the nonlinear solvers.
789 
790    If the user did did not write their equations in implicit form, this
791    function recasts them in implicit form.
792 
793    Level: developer
794 
795 .keywords: TS, compute
796 
797 .seealso: TSSetIFunction(), TSComputeRHSFunction()
798 @*/
799 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec Y,PetscBool imex)
800 {
801   PetscErrorCode ierr;
802   TSIFunction    ifunction;
803   TSRHSFunction  rhsfunction;
804   void           *ctx;
805   DM             dm;
806 
807   PetscFunctionBegin;
808   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
809   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
810   PetscValidHeaderSpecific(Udot,VEC_CLASSID,4);
811   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
812 
813   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
814   ierr = DMTSGetIFunction(dm,&ifunction,&ctx);CHKERRQ(ierr);
815   ierr = DMTSGetRHSFunction(dm,&rhsfunction,NULL);CHKERRQ(ierr);
816 
817   if (!rhsfunction && !ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
818 
819   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Udot,Y);CHKERRQ(ierr);
820   if (ifunction) {
821     PetscStackPush("TS user implicit function");
822     ierr = (*ifunction)(ts,t,U,Udot,Y,ctx);CHKERRQ(ierr);
823     PetscStackPop;
824   }
825   if (imex) {
826     if (!ifunction) {
827       ierr = VecCopy(Udot,Y);CHKERRQ(ierr);
828     }
829   } else if (rhsfunction) {
830     if (ifunction) {
831       Vec Frhs;
832       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
833       ierr = TSComputeRHSFunction(ts,t,U,Frhs);CHKERRQ(ierr);
834       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
835     } else {
836       ierr = TSComputeRHSFunction(ts,t,U,Y);CHKERRQ(ierr);
837       ierr = VecAYPX(Y,-1,Udot);CHKERRQ(ierr);
838     }
839   }
840   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Udot,Y);CHKERRQ(ierr);
841   PetscFunctionReturn(0);
842 }
843 
844 #undef __FUNCT__
845 #define __FUNCT__ "TSComputeIJacobian"
846 /*@
847    TSComputeIJacobian - Evaluates the Jacobian of the DAE
848 
849    Collective on TS and Vec
850 
851    Input
852       Input Parameters:
853 +  ts - the TS context
854 .  t - current timestep
855 .  U - state vector
856 .  Udot - time derivative of state vector
857 .  shift - shift to apply, see note below
858 -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
859 
860    Output Parameters:
861 +  A - Jacobian matrix
862 .  B - optional preconditioning matrix
863 -  flag - flag indicating matrix structure
864 
865    Notes:
866    If F(t,U,Udot)=0 is the DAE, the required Jacobian is
867 
868    dF/dU + shift*dF/dUdot
869 
870    Most users should not need to explicitly call this routine, as it
871    is used internally within the nonlinear solvers.
872 
873    Level: developer
874 
875 .keywords: TS, compute, Jacobian, matrix
876 
877 .seealso:  TSSetIJacobian()
878 @*/
879 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,PetscBool imex)
880 {
881   PetscErrorCode ierr;
882   TSIJacobian    ijacobian;
883   TSRHSJacobian  rhsjacobian;
884   DM             dm;
885   void           *ctx;
886 
887   PetscFunctionBegin;
888   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
889   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
890   PetscValidHeaderSpecific(Udot,VEC_CLASSID,4);
891   PetscValidPointer(A,6);
892   PetscValidHeaderSpecific(A,MAT_CLASSID,6);
893   PetscValidPointer(B,7);
894   PetscValidHeaderSpecific(B,MAT_CLASSID,7);
895 
896   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
897   ierr = DMTSGetIJacobian(dm,&ijacobian,&ctx);CHKERRQ(ierr);
898   ierr = DMTSGetRHSJacobian(dm,&rhsjacobian,NULL);CHKERRQ(ierr);
899 
900   if (!rhsjacobian && !ijacobian) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
901 
902   ierr = PetscLogEventBegin(TS_JacobianEval,ts,U,A,B);CHKERRQ(ierr);
903   if (ijacobian) {
904     PetscBool missing;
905     PetscStackPush("TS user implicit Jacobian");
906     ierr = (*ijacobian)(ts,t,U,Udot,shift,A,B,ctx);CHKERRQ(ierr);
907     PetscStackPop;
908     if (A) {
909       ierr = MatMissingDiagonal(A,&missing,NULL);CHKERRQ(ierr);
910       if (missing) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Amat passed to TSSetIJacobian() must have all diagonal entries set, if they are zero you must still set them with a zero value");
911     }
912     if (B && B != A) {
913       ierr = MatMissingDiagonal(B,&missing,NULL);CHKERRQ(ierr);
914       if (missing) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Bmat passed to TSSetIJacobian() must have all diagonal entries set, if they are zero you must still set them with a zero value");
915     }
916   }
917   if (imex) {
918     if (!ijacobian) {  /* system was written as Udot = G(t,U) */
919       ierr = MatZeroEntries(A);CHKERRQ(ierr);
920       ierr = MatShift(A,shift);CHKERRQ(ierr);
921       if (A != B) {
922         ierr = MatZeroEntries(B);CHKERRQ(ierr);
923         ierr = MatShift(B,shift);CHKERRQ(ierr);
924       }
925     }
926   } else {
927     Mat Arhs = NULL,Brhs = NULL;
928     if (rhsjacobian) {
929       if (ijacobian) {
930         ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
931       } else {
932         ierr = TSGetIJacobian(ts,&Arhs,&Brhs,NULL,NULL);CHKERRQ(ierr);
933       }
934       ierr = TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);CHKERRQ(ierr);
935     }
936     if (Arhs == A) {           /* No IJacobian, so we only have the RHS matrix */
937       ts->rhsjacobian.scale = -1;
938       ts->rhsjacobian.shift = shift;
939       ierr = MatScale(A,-1);CHKERRQ(ierr);
940       ierr = MatShift(A,shift);CHKERRQ(ierr);
941       if (A != B) {
942         ierr = MatScale(B,-1);CHKERRQ(ierr);
943         ierr = MatShift(B,shift);CHKERRQ(ierr);
944       }
945     } else if (Arhs) {          /* Both IJacobian and RHSJacobian */
946       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
947       if (!ijacobian) {         /* No IJacobian provided, but we have a separate RHS matrix */
948         ierr = MatZeroEntries(A);CHKERRQ(ierr);
949         ierr = MatShift(A,shift);CHKERRQ(ierr);
950         if (A != B) {
951           ierr = MatZeroEntries(B);CHKERRQ(ierr);
952           ierr = MatShift(B,shift);CHKERRQ(ierr);
953         }
954       }
955       ierr = MatAXPY(A,-1,Arhs,axpy);CHKERRQ(ierr);
956       if (A != B) {
957         ierr = MatAXPY(B,-1,Brhs,axpy);CHKERRQ(ierr);
958       }
959     }
960   }
961   ierr = PetscLogEventEnd(TS_JacobianEval,ts,U,A,B);CHKERRQ(ierr);
962   PetscFunctionReturn(0);
963 }
964 
965 #undef __FUNCT__
966 #define __FUNCT__ "TSSetRHSFunction"
967 /*@C
968     TSSetRHSFunction - Sets the routine for evaluating the function,
969     where U_t = G(t,u).
970 
971     Logically Collective on TS
972 
973     Input Parameters:
974 +   ts - the TS context obtained from TSCreate()
975 .   r - vector to put the computed right hand side (or NULL to have it created)
976 .   f - routine for evaluating the right-hand-side function
977 -   ctx - [optional] user-defined context for private data for the
978           function evaluation routine (may be NULL)
979 
980     Calling sequence of func:
981 $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
982 
983 +   t - current timestep
984 .   u - input vector
985 .   F - function vector
986 -   ctx - [optional] user-defined function context
987 
988     Level: beginner
989 
990     Notes: You must call this function or TSSetIFunction() to define your ODE. You cannot use this function when solving a DAE.
991 
992 .keywords: TS, timestep, set, right-hand-side, function
993 
994 .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSSetIFunction()
995 @*/
996 PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
997 {
998   PetscErrorCode ierr;
999   SNES           snes;
1000   Vec            ralloc = NULL;
1001   DM             dm;
1002 
1003   PetscFunctionBegin;
1004   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1005   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
1006 
1007   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1008   ierr = DMTSSetRHSFunction(dm,f,ctx);CHKERRQ(ierr);
1009   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1010   if (!r && !ts->dm && ts->vec_sol) {
1011     ierr = VecDuplicate(ts->vec_sol,&ralloc);CHKERRQ(ierr);
1012     r    = ralloc;
1013   }
1014   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
1015   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
1016   PetscFunctionReturn(0);
1017 }
1018 
1019 #undef __FUNCT__
1020 #define __FUNCT__ "TSSetSolutionFunction"
1021 /*@C
1022     TSSetSolutionFunction - Provide a function that computes the solution of the ODE or DAE
1023 
1024     Logically Collective on TS
1025 
1026     Input Parameters:
1027 +   ts - the TS context obtained from TSCreate()
1028 .   f - routine for evaluating the solution
1029 -   ctx - [optional] user-defined context for private data for the
1030           function evaluation routine (may be NULL)
1031 
1032     Calling sequence of func:
1033 $     func (TS ts,PetscReal t,Vec u,void *ctx);
1034 
1035 +   t - current timestep
1036 .   u - output vector
1037 -   ctx - [optional] user-defined function context
1038 
1039     Notes:
1040     This routine is used for testing accuracy of time integration schemes when you already know the solution.
1041     If analytic solutions are not known for your system, consider using the Method of Manufactured Solutions to
1042     create closed-form solutions with non-physical forcing terms.
1043 
1044     For low-dimensional problems solved in serial, such as small discrete systems, TSMonitorLGError() can be used to monitor the error history.
1045 
1046     Level: beginner
1047 
1048 .keywords: TS, timestep, set, right-hand-side, function
1049 
1050 .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetForcingFunction()
1051 @*/
1052 PetscErrorCode  TSSetSolutionFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
1053 {
1054   PetscErrorCode ierr;
1055   DM             dm;
1056 
1057   PetscFunctionBegin;
1058   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1059   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1060   ierr = DMTSSetSolutionFunction(dm,f,ctx);CHKERRQ(ierr);
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 #undef __FUNCT__
1065 #define __FUNCT__ "TSSetForcingFunction"
1066 /*@C
1067     TSSetForcingFunction - Provide a function that computes a forcing term for a ODE or PDE
1068 
1069     Logically Collective on TS
1070 
1071     Input Parameters:
1072 +   ts - the TS context obtained from TSCreate()
1073 .   f - routine for evaluating the forcing function
1074 -   ctx - [optional] user-defined context for private data for the
1075           function evaluation routine (may be NULL)
1076 
1077     Calling sequence of func:
1078 $     func (TS ts,PetscReal t,Vec u,void *ctx);
1079 
1080 +   t - current timestep
1081 .   u - output vector
1082 -   ctx - [optional] user-defined function context
1083 
1084     Notes:
1085     This routine is useful for testing accuracy of time integration schemes when using the Method of Manufactured Solutions to
1086     create closed-form solutions with a non-physical forcing term.
1087 
1088     For low-dimensional problems solved in serial, such as small discrete systems, TSMonitorLGError() can be used to monitor the error history.
1089 
1090     Level: beginner
1091 
1092 .keywords: TS, timestep, set, right-hand-side, function
1093 
1094 .seealso: TSSetRHSJacobian(), TSSetIJacobian(), TSComputeSolutionFunction(), TSSetSolutionFunction()
1095 @*/
1096 PetscErrorCode  TSSetForcingFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
1097 {
1098   PetscErrorCode ierr;
1099   DM             dm;
1100 
1101   PetscFunctionBegin;
1102   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1103   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1104   ierr = DMTSSetForcingFunction(dm,f,ctx);CHKERRQ(ierr);
1105   PetscFunctionReturn(0);
1106 }
1107 
1108 #undef __FUNCT__
1109 #define __FUNCT__ "TSSetRHSJacobian"
1110 /*@C
1111    TSSetRHSJacobian - Sets the function to compute the Jacobian of G,
1112    where U_t = G(U,t), as well as the location to store the matrix.
1113 
1114    Logically Collective on TS
1115 
1116    Input Parameters:
1117 +  ts  - the TS context obtained from TSCreate()
1118 .  Amat - (approximate) Jacobian matrix
1119 .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1120 .  f   - the Jacobian evaluation routine
1121 -  ctx - [optional] user-defined context for private data for the
1122          Jacobian evaluation routine (may be NULL)
1123 
1124    Calling sequence of f:
1125 $     func (TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx);
1126 
1127 +  t - current timestep
1128 .  u - input vector
1129 .  Amat - (approximate) Jacobian matrix
1130 .  Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat)
1131 -  ctx - [optional] user-defined context for matrix evaluation routine
1132 
1133    Notes:
1134    You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value
1135 
1136    The TS solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f()
1137    You should not assume the values are the same in the next call to f() as you set them in the previous call.
1138 
1139    Level: beginner
1140 
1141 .keywords: TS, timestep, set, right-hand-side, Jacobian
1142 
1143 .seealso: SNESComputeJacobianDefaultColor(), TSSetRHSFunction(), TSRHSJacobianSetReuse(), TSSetIJacobian()
1144 
1145 @*/
1146 PetscErrorCode  TSSetRHSJacobian(TS ts,Mat Amat,Mat Pmat,TSRHSJacobian f,void *ctx)
1147 {
1148   PetscErrorCode ierr;
1149   SNES           snes;
1150   DM             dm;
1151   TSIJacobian    ijacobian;
1152 
1153   PetscFunctionBegin;
1154   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1155   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1156   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
1157   if (Amat) PetscCheckSameComm(ts,1,Amat,2);
1158   if (Pmat) PetscCheckSameComm(ts,1,Pmat,3);
1159 
1160   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1161   ierr = DMTSSetRHSJacobian(dm,f,ctx);CHKERRQ(ierr);
1162   if (f == TSComputeRHSJacobianConstant) {
1163     /* Handle this case automatically for the user; otherwise user should call themselves. */
1164     ierr = TSRHSJacobianSetReuse(ts,PETSC_TRUE);CHKERRQ(ierr);
1165   }
1166   ierr = DMTSGetIJacobian(dm,&ijacobian,NULL);CHKERRQ(ierr);
1167   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1168   if (!ijacobian) {
1169     ierr = SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);CHKERRQ(ierr);
1170   }
1171   if (Amat) {
1172     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
1173     ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
1174 
1175     ts->Arhs = Amat;
1176   }
1177   if (Pmat) {
1178     ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);
1179     ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1180 
1181     ts->Brhs = Pmat;
1182   }
1183   PetscFunctionReturn(0);
1184 }
1185 
1186 
1187 #undef __FUNCT__
1188 #define __FUNCT__ "TSSetIFunction"
1189 /*@C
1190    TSSetIFunction - Set the function to compute F(t,U,U_t) where F() = 0 is the DAE to be solved.
1191 
1192    Logically Collective on TS
1193 
1194    Input Parameters:
1195 +  ts  - the TS context obtained from TSCreate()
1196 .  r   - vector to hold the residual (or NULL to have it created internally)
1197 .  f   - the function evaluation routine
1198 -  ctx - user-defined context for private data for the function evaluation routine (may be NULL)
1199 
1200    Calling sequence of f:
1201 $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
1202 
1203 +  t   - time at step/stage being solved
1204 .  u   - state vector
1205 .  u_t - time derivative of state vector
1206 .  F   - function vector
1207 -  ctx - [optional] user-defined context for matrix evaluation routine
1208 
1209    Important:
1210    The user MUST call either this routine or TSSetRHSFunction() to define the ODE.  When solving DAEs you must use this function.
1211 
1212    Level: beginner
1213 
1214 .keywords: TS, timestep, set, DAE, Jacobian
1215 
1216 .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
1217 @*/
1218 PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
1219 {
1220   PetscErrorCode ierr;
1221   SNES           snes;
1222   Vec            resalloc = NULL;
1223   DM             dm;
1224 
1225   PetscFunctionBegin;
1226   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1227   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
1228 
1229   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1230   ierr = DMTSSetIFunction(dm,f,ctx);CHKERRQ(ierr);
1231 
1232   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1233   if (!res && !ts->dm && ts->vec_sol) {
1234     ierr = VecDuplicate(ts->vec_sol,&resalloc);CHKERRQ(ierr);
1235     res  = resalloc;
1236   }
1237   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
1238   ierr = VecDestroy(&resalloc);CHKERRQ(ierr);
1239   PetscFunctionReturn(0);
1240 }
1241 
1242 #undef __FUNCT__
1243 #define __FUNCT__ "TSGetIFunction"
1244 /*@C
1245    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
1246 
1247    Not Collective
1248 
1249    Input Parameter:
1250 .  ts - the TS context
1251 
1252    Output Parameter:
1253 +  r - vector to hold residual (or NULL)
1254 .  func - the function to compute residual (or NULL)
1255 -  ctx - the function context (or NULL)
1256 
1257    Level: advanced
1258 
1259 .keywords: TS, nonlinear, get, function
1260 
1261 .seealso: TSSetIFunction(), SNESGetFunction()
1262 @*/
1263 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
1264 {
1265   PetscErrorCode ierr;
1266   SNES           snes;
1267   DM             dm;
1268 
1269   PetscFunctionBegin;
1270   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1271   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1272   ierr = SNESGetFunction(snes,r,NULL,NULL);CHKERRQ(ierr);
1273   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1274   ierr = DMTSGetIFunction(dm,func,ctx);CHKERRQ(ierr);
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 #undef __FUNCT__
1279 #define __FUNCT__ "TSGetRHSFunction"
1280 /*@C
1281    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
1282 
1283    Not Collective
1284 
1285    Input Parameter:
1286 .  ts - the TS context
1287 
1288    Output Parameter:
1289 +  r - vector to hold computed right hand side (or NULL)
1290 .  func - the function to compute right hand side (or NULL)
1291 -  ctx - the function context (or NULL)
1292 
1293    Level: advanced
1294 
1295 .keywords: TS, nonlinear, get, function
1296 
1297 .seealso: TSSetRHSFunction(), SNESGetFunction()
1298 @*/
1299 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
1300 {
1301   PetscErrorCode ierr;
1302   SNES           snes;
1303   DM             dm;
1304 
1305   PetscFunctionBegin;
1306   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1307   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1308   ierr = SNESGetFunction(snes,r,NULL,NULL);CHKERRQ(ierr);
1309   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1310   ierr = DMTSGetRHSFunction(dm,func,ctx);CHKERRQ(ierr);
1311   PetscFunctionReturn(0);
1312 }
1313 
1314 #undef __FUNCT__
1315 #define __FUNCT__ "TSSetIJacobian"
1316 /*@C
1317    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
1318         provided with TSSetIFunction().
1319 
1320    Logically Collective on TS
1321 
1322    Input Parameters:
1323 +  ts  - the TS context obtained from TSCreate()
1324 .  Amat - (approximate) Jacobian matrix
1325 .  Pmat - matrix used to compute preconditioner (usually the same as Amat)
1326 .  f   - the Jacobian evaluation routine
1327 -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be NULL)
1328 
1329    Calling sequence of f:
1330 $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx);
1331 
1332 +  t    - time at step/stage being solved
1333 .  U    - state vector
1334 .  U_t  - time derivative of state vector
1335 .  a    - shift
1336 .  Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
1337 .  Pmat - matrix used for constructing preconditioner, usually the same as Amat
1338 -  ctx  - [optional] user-defined context for matrix evaluation routine
1339 
1340    Notes:
1341    The matrices Amat and Pmat are exactly the matrices that are used by SNES for the nonlinear solve.
1342 
1343    If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
1344    space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
1345 
1346    The matrix dF/dU + a*dF/dU_t you provide turns out to be
1347    the Jacobian of F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
1348    The time integrator internally approximates U_t by W+a*U where the positive "shift"
1349    a and vector W depend on the integration method, step size, and past states. For example with
1350    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
1351    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
1352 
1353    You must set all the diagonal entries of the matrices, if they are zero you must still set them with a zero value
1354 
1355    The TS solver may modify the nonzero structure and the entries of the matrices Amat and Pmat between the calls to f()
1356    You should not assume the values are the same in the next call to f() as you set them in the previous call.
1357 
1358    Level: beginner
1359 
1360 .keywords: TS, timestep, DAE, Jacobian
1361 
1362 .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESComputeJacobianDefaultColor(), SNESComputeJacobianDefault(), TSSetRHSFunction()
1363 
1364 @*/
1365 PetscErrorCode  TSSetIJacobian(TS ts,Mat Amat,Mat Pmat,TSIJacobian f,void *ctx)
1366 {
1367   PetscErrorCode ierr;
1368   SNES           snes;
1369   DM             dm;
1370 
1371   PetscFunctionBegin;
1372   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1373   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1374   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
1375   if (Amat) PetscCheckSameComm(ts,1,Amat,2);
1376   if (Pmat) PetscCheckSameComm(ts,1,Pmat,3);
1377 
1378   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1379   ierr = DMTSSetIJacobian(dm,f,ctx);CHKERRQ(ierr);
1380 
1381   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1382   ierr = SNESSetJacobian(snes,Amat,Pmat,SNESTSFormJacobian,ts);CHKERRQ(ierr);
1383   PetscFunctionReturn(0);
1384 }
1385 
1386 #undef __FUNCT__
1387 #define __FUNCT__ "TSRHSJacobianSetReuse"
1388 /*@
1389    TSRHSJacobianSetReuse - restore RHS Jacobian before re-evaluating.  Without this flag, TS will change the sign and
1390    shift the RHS Jacobian for a finite-time-step implicit solve, in which case the user function will need to recompute
1391    the entire Jacobian.  The reuse flag must be set if the evaluation function will assume that the matrix entries have
1392    not been changed by the TS.
1393 
1394    Logically Collective
1395 
1396    Input Arguments:
1397 +  ts - TS context obtained from TSCreate()
1398 -  reuse - PETSC_TRUE if the RHS Jacobian
1399 
1400    Level: intermediate
1401 
1402 .seealso: TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
1403 @*/
1404 PetscErrorCode TSRHSJacobianSetReuse(TS ts,PetscBool reuse)
1405 {
1406   PetscFunctionBegin;
1407   ts->rhsjacobian.reuse = reuse;
1408   PetscFunctionReturn(0);
1409 }
1410 
1411 #undef __FUNCT__
1412 #define __FUNCT__ "TSLoad"
1413 /*@C
1414   TSLoad - Loads a KSP that has been stored in binary  with KSPView().
1415 
1416   Collective on PetscViewer
1417 
1418   Input Parameters:
1419 + newdm - the newly loaded TS, this needs to have been created with TSCreate() or
1420            some related function before a call to TSLoad().
1421 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
1422 
1423    Level: intermediate
1424 
1425   Notes:
1426    The type is determined by the data in the file, any type set into the TS before this call is ignored.
1427 
1428   Notes for advanced users:
1429   Most users should not need to know the details of the binary storage
1430   format, since TSLoad() and TSView() completely hide these details.
1431   But for anyone who's interested, the standard binary matrix storage
1432   format is
1433 .vb
1434      has not yet been determined
1435 .ve
1436 
1437 .seealso: PetscViewerBinaryOpen(), TSView(), MatLoad(), VecLoad()
1438 @*/
1439 PetscErrorCode  TSLoad(TS ts, PetscViewer viewer)
1440 {
1441   PetscErrorCode ierr;
1442   PetscBool      isbinary;
1443   PetscInt       classid;
1444   char           type[256];
1445   DMTS           sdm;
1446   DM             dm;
1447 
1448   PetscFunctionBegin;
1449   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1450   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1451   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1452   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1453 
1454   ierr = PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);CHKERRQ(ierr);
1455   if (classid != TS_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONG,"Not TS next in file");
1456   ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr);
1457   ierr = TSSetType(ts, type);CHKERRQ(ierr);
1458   if (ts->ops->load) {
1459     ierr = (*ts->ops->load)(ts,viewer);CHKERRQ(ierr);
1460   }
1461   ierr = DMCreate(PetscObjectComm((PetscObject)ts),&dm);CHKERRQ(ierr);
1462   ierr = DMLoad(dm,viewer);CHKERRQ(ierr);
1463   ierr = TSSetDM(ts,dm);CHKERRQ(ierr);
1464   ierr = DMCreateGlobalVector(ts->dm,&ts->vec_sol);CHKERRQ(ierr);
1465   ierr = VecLoad(ts->vec_sol,viewer);CHKERRQ(ierr);
1466   ierr = DMGetDMTS(ts->dm,&sdm);CHKERRQ(ierr);
1467   ierr = DMTSLoad(sdm,viewer);CHKERRQ(ierr);
1468   PetscFunctionReturn(0);
1469 }
1470 
1471 #include <petscdraw.h>
1472 #if defined(PETSC_HAVE_SAWS)
1473 #include <petscviewersaws.h>
1474 #endif
1475 #undef __FUNCT__
1476 #define __FUNCT__ "TSView"
1477 /*@C
1478     TSView - Prints the TS data structure.
1479 
1480     Collective on TS
1481 
1482     Input Parameters:
1483 +   ts - the TS context obtained from TSCreate()
1484 -   viewer - visualization context
1485 
1486     Options Database Key:
1487 .   -ts_view - calls TSView() at end of TSStep()
1488 
1489     Notes:
1490     The available visualization contexts include
1491 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1492 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1493          output where only the first processor opens
1494          the file.  All other processors send their
1495          data to the first processor to print.
1496 
1497     The user can open an alternative visualization context with
1498     PetscViewerASCIIOpen() - output to a specified file.
1499 
1500     Level: beginner
1501 
1502 .keywords: TS, timestep, view
1503 
1504 .seealso: PetscViewerASCIIOpen()
1505 @*/
1506 PetscErrorCode  TSView(TS ts,PetscViewer viewer)
1507 {
1508   PetscErrorCode ierr;
1509   TSType         type;
1510   PetscBool      iascii,isstring,isundials,isbinary,isdraw;
1511   DMTS           sdm;
1512 #if defined(PETSC_HAVE_SAWS)
1513   PetscBool      issaws;
1514 #endif
1515 
1516   PetscFunctionBegin;
1517   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1518   if (!viewer) {
1519     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ts),&viewer);CHKERRQ(ierr);
1520   }
1521   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1522   PetscCheckSameComm(ts,1,viewer,2);
1523 
1524   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1525   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
1526   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
1527   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1528 #if defined(PETSC_HAVE_SAWS)
1529   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);CHKERRQ(ierr);
1530 #endif
1531   if (iascii) {
1532     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer);CHKERRQ(ierr);
1533     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
1534     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",(double)ts->max_time);CHKERRQ(ierr);
1535     if (ts->problem_type == TS_NONLINEAR) {
1536       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);CHKERRQ(ierr);
1537       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr);
1538     }
1539     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);CHKERRQ(ierr);
1540     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
1541     ierr = DMGetDMTS(ts->dm,&sdm);CHKERRQ(ierr);
1542     ierr = DMTSView(sdm,viewer);CHKERRQ(ierr);
1543     if (ts->ops->view) {
1544       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1545       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1546       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1547     }
1548   } else if (isstring) {
1549     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
1550     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
1551   } else if (isbinary) {
1552     PetscInt    classid = TS_FILE_CLASSID;
1553     MPI_Comm    comm;
1554     PetscMPIInt rank;
1555     char        type[256];
1556 
1557     ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1558     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1559     if (!rank) {
1560       ierr = PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
1561       ierr = PetscStrncpy(type,((PetscObject)ts)->type_name,256);CHKERRQ(ierr);
1562       ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
1563     }
1564     if (ts->ops->view) {
1565       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1566     }
1567     ierr = DMView(ts->dm,viewer);CHKERRQ(ierr);
1568     ierr = VecView(ts->vec_sol,viewer);CHKERRQ(ierr);
1569     ierr = DMGetDMTS(ts->dm,&sdm);CHKERRQ(ierr);
1570     ierr = DMTSView(sdm,viewer);CHKERRQ(ierr);
1571   } else if (isdraw) {
1572     PetscDraw draw;
1573     char      str[36];
1574     PetscReal x,y,bottom,h;
1575 
1576     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
1577     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
1578     ierr   = PetscStrcpy(str,"TS: ");CHKERRQ(ierr);
1579     ierr   = PetscStrcat(str,((PetscObject)ts)->type_name);CHKERRQ(ierr);
1580     ierr   = PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,str,NULL,&h);CHKERRQ(ierr);
1581     bottom = y - h;
1582     ierr   = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
1583     if (ts->ops->view) {
1584       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1585     }
1586     ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
1587 #if defined(PETSC_HAVE_SAWS)
1588   } else if (issaws) {
1589     PetscMPIInt rank;
1590     const char  *name;
1591 
1592     ierr = PetscObjectGetName((PetscObject)ts,&name);CHKERRQ(ierr);
1593     ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
1594     if (!((PetscObject)ts)->amsmem && !rank) {
1595       char       dir[1024];
1596 
1597       ierr = PetscObjectViewSAWs((PetscObject)ts,viewer);CHKERRQ(ierr);
1598       ierr = PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time_step",name);CHKERRQ(ierr);
1599       PetscStackCallSAWs(SAWs_Register,(dir,&ts->steps,1,SAWs_READ,SAWs_INT));
1600       ierr = PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/time",name);CHKERRQ(ierr);
1601       PetscStackCallSAWs(SAWs_Register,(dir,&ts->ptime,1,SAWs_READ,SAWs_DOUBLE));
1602     }
1603     if (ts->ops->view) {
1604       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1605     }
1606 #endif
1607   }
1608 
1609   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1610   ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
1611   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1612   PetscFunctionReturn(0);
1613 }
1614 
1615 
1616 #undef __FUNCT__
1617 #define __FUNCT__ "TSSetApplicationContext"
1618 /*@
1619    TSSetApplicationContext - Sets an optional user-defined context for
1620    the timesteppers.
1621 
1622    Logically Collective on TS
1623 
1624    Input Parameters:
1625 +  ts - the TS context obtained from TSCreate()
1626 -  usrP - optional user context
1627 
1628    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
1629     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1630 
1631    Level: intermediate
1632 
1633 .keywords: TS, timestep, set, application, context
1634 
1635 .seealso: TSGetApplicationContext()
1636 @*/
1637 PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
1638 {
1639   PetscFunctionBegin;
1640   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1641   ts->user = usrP;
1642   PetscFunctionReturn(0);
1643 }
1644 
1645 #undef __FUNCT__
1646 #define __FUNCT__ "TSGetApplicationContext"
1647 /*@
1648     TSGetApplicationContext - Gets the user-defined context for the
1649     timestepper.
1650 
1651     Not Collective
1652 
1653     Input Parameter:
1654 .   ts - the TS context obtained from TSCreate()
1655 
1656     Output Parameter:
1657 .   usrP - user context
1658 
1659    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
1660     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1661 
1662     Level: intermediate
1663 
1664 .keywords: TS, timestep, get, application, context
1665 
1666 .seealso: TSSetApplicationContext()
1667 @*/
1668 PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
1669 {
1670   PetscFunctionBegin;
1671   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1672   *(void**)usrP = ts->user;
1673   PetscFunctionReturn(0);
1674 }
1675 
1676 #undef __FUNCT__
1677 #define __FUNCT__ "TSGetTimeStepNumber"
1678 /*@
1679    TSGetTimeStepNumber - Gets the number of time steps completed.
1680 
1681    Not Collective
1682 
1683    Input Parameter:
1684 .  ts - the TS context obtained from TSCreate()
1685 
1686    Output Parameter:
1687 .  iter - number of steps completed so far
1688 
1689    Level: intermediate
1690 
1691 .keywords: TS, timestep, get, iteration, number
1692 .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSSetPostStep()
1693 @*/
1694 PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt *iter)
1695 {
1696   PetscFunctionBegin;
1697   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1698   PetscValidIntPointer(iter,2);
1699   *iter = ts->steps;
1700   PetscFunctionReturn(0);
1701 }
1702 
1703 #undef __FUNCT__
1704 #define __FUNCT__ "TSSetInitialTimeStep"
1705 /*@
1706    TSSetInitialTimeStep - Sets the initial timestep to be used,
1707    as well as the initial time.
1708 
1709    Logically Collective on TS
1710 
1711    Input Parameters:
1712 +  ts - the TS context obtained from TSCreate()
1713 .  initial_time - the initial time
1714 -  time_step - the size of the timestep
1715 
1716    Level: intermediate
1717 
1718 .seealso: TSSetTimeStep(), TSGetTimeStep()
1719 
1720 .keywords: TS, set, initial, timestep
1721 @*/
1722 PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
1723 {
1724   PetscErrorCode ierr;
1725 
1726   PetscFunctionBegin;
1727   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1728   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
1729   ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr);
1730   PetscFunctionReturn(0);
1731 }
1732 
1733 #undef __FUNCT__
1734 #define __FUNCT__ "TSSetTimeStep"
1735 /*@
1736    TSSetTimeStep - Allows one to reset the timestep at any time,
1737    useful for simple pseudo-timestepping codes.
1738 
1739    Logically Collective on TS
1740 
1741    Input Parameters:
1742 +  ts - the TS context obtained from TSCreate()
1743 -  time_step - the size of the timestep
1744 
1745    Level: intermediate
1746 
1747 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1748 
1749 .keywords: TS, set, timestep
1750 @*/
1751 PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
1752 {
1753   PetscFunctionBegin;
1754   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1755   PetscValidLogicalCollectiveReal(ts,time_step,2);
1756   ts->time_step      = time_step;
1757   ts->time_step_orig = time_step;
1758   PetscFunctionReturn(0);
1759 }
1760 
1761 #undef __FUNCT__
1762 #define __FUNCT__ "TSSetExactFinalTime"
1763 /*@
1764    TSSetExactFinalTime - Determines whether to adapt the final time step to
1765      match the exact final time, interpolate solution to the exact final time,
1766      or just return at the final time TS computed.
1767 
1768   Logically Collective on TS
1769 
1770    Input Parameter:
1771 +   ts - the time-step context
1772 -   eftopt - exact final time option
1773 
1774    Level: beginner
1775 
1776 .seealso: TSExactFinalTimeOption
1777 @*/
1778 PetscErrorCode  TSSetExactFinalTime(TS ts,TSExactFinalTimeOption eftopt)
1779 {
1780   PetscFunctionBegin;
1781   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1782   PetscValidLogicalCollectiveEnum(ts,eftopt,2);
1783   ts->exact_final_time = eftopt;
1784   PetscFunctionReturn(0);
1785 }
1786 
1787 #undef __FUNCT__
1788 #define __FUNCT__ "TSGetTimeStep"
1789 /*@
1790    TSGetTimeStep - Gets the current timestep size.
1791 
1792    Not Collective
1793 
1794    Input Parameter:
1795 .  ts - the TS context obtained from TSCreate()
1796 
1797    Output Parameter:
1798 .  dt - the current timestep size
1799 
1800    Level: intermediate
1801 
1802 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1803 
1804 .keywords: TS, get, timestep
1805 @*/
1806 PetscErrorCode  TSGetTimeStep(TS ts,PetscReal *dt)
1807 {
1808   PetscFunctionBegin;
1809   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1810   PetscValidRealPointer(dt,2);
1811   *dt = ts->time_step;
1812   PetscFunctionReturn(0);
1813 }
1814 
1815 #undef __FUNCT__
1816 #define __FUNCT__ "TSGetSolution"
1817 /*@
1818    TSGetSolution - Returns the solution at the present timestep. It
1819    is valid to call this routine inside the function that you are evaluating
1820    in order to move to the new timestep. This vector not changed until
1821    the solution at the next timestep has been calculated.
1822 
1823    Not Collective, but Vec returned is parallel if TS is parallel
1824 
1825    Input Parameter:
1826 .  ts - the TS context obtained from TSCreate()
1827 
1828    Output Parameter:
1829 .  v - the vector containing the solution
1830 
1831    Level: intermediate
1832 
1833 .seealso: TSGetTimeStep()
1834 
1835 .keywords: TS, timestep, get, solution
1836 @*/
1837 PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1838 {
1839   PetscFunctionBegin;
1840   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1841   PetscValidPointer(v,2);
1842   *v = ts->vec_sol;
1843   PetscFunctionReturn(0);
1844 }
1845 
1846 #undef __FUNCT__
1847 #define __FUNCT__ "TSGetCostGradients"
1848 /*@
1849    TSGetCostGradients - Returns the gradients from the TSAdjointSolve()
1850 
1851    Not Collective, but Vec returned is parallel if TS is parallel
1852 
1853    Input Parameter:
1854 .  ts - the TS context obtained from TSCreate()
1855 
1856    Output Parameter:
1857 +  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1858 -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
1859 
1860    Level: intermediate
1861 
1862 .seealso: TSGetTimeStep()
1863 
1864 .keywords: TS, timestep, get, sensitivity
1865 @*/
1866 PetscErrorCode  TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu)
1867 {
1868   PetscFunctionBegin;
1869   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1870   if (numcost) *numcost = ts->numcost;
1871   if (lambda)  *lambda  = ts->vecs_sensi;
1872   if (mu)      *mu      = ts->vecs_sensip;
1873   PetscFunctionReturn(0);
1874 }
1875 
1876 /* ----- Routines to initialize and destroy a timestepper ---- */
1877 #undef __FUNCT__
1878 #define __FUNCT__ "TSSetProblemType"
1879 /*@
1880   TSSetProblemType - Sets the type of problem to be solved.
1881 
1882   Not collective
1883 
1884   Input Parameters:
1885 + ts   - The TS
1886 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1887 .vb
1888          U_t - A U = 0      (linear)
1889          U_t - A(t) U = 0   (linear)
1890          F(t,U,U_t) = 0     (nonlinear)
1891 .ve
1892 
1893    Level: beginner
1894 
1895 .keywords: TS, problem type
1896 .seealso: TSSetUp(), TSProblemType, TS
1897 @*/
1898 PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1899 {
1900   PetscErrorCode ierr;
1901 
1902   PetscFunctionBegin;
1903   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1904   ts->problem_type = type;
1905   if (type == TS_LINEAR) {
1906     SNES snes;
1907     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1908     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
1909   }
1910   PetscFunctionReturn(0);
1911 }
1912 
1913 #undef __FUNCT__
1914 #define __FUNCT__ "TSGetProblemType"
1915 /*@C
1916   TSGetProblemType - Gets the type of problem to be solved.
1917 
1918   Not collective
1919 
1920   Input Parameter:
1921 . ts   - The TS
1922 
1923   Output Parameter:
1924 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1925 .vb
1926          M U_t = A U
1927          M(t) U_t = A(t) U
1928          F(t,U,U_t)
1929 .ve
1930 
1931    Level: beginner
1932 
1933 .keywords: TS, problem type
1934 .seealso: TSSetUp(), TSProblemType, TS
1935 @*/
1936 PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1937 {
1938   PetscFunctionBegin;
1939   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1940   PetscValidIntPointer(type,2);
1941   *type = ts->problem_type;
1942   PetscFunctionReturn(0);
1943 }
1944 
1945 #undef __FUNCT__
1946 #define __FUNCT__ "TSSetUp"
1947 /*@
1948    TSSetUp - Sets up the internal data structures for the later use
1949    of a timestepper.
1950 
1951    Collective on TS
1952 
1953    Input Parameter:
1954 .  ts - the TS context obtained from TSCreate()
1955 
1956    Notes:
1957    For basic use of the TS solvers the user need not explicitly call
1958    TSSetUp(), since these actions will automatically occur during
1959    the call to TSStep().  However, if one wishes to control this
1960    phase separately, TSSetUp() should be called after TSCreate()
1961    and optional routines of the form TSSetXXX(), but before TSStep().
1962 
1963    Level: advanced
1964 
1965 .keywords: TS, timestep, setup
1966 
1967 .seealso: TSCreate(), TSStep(), TSDestroy()
1968 @*/
1969 PetscErrorCode  TSSetUp(TS ts)
1970 {
1971   PetscErrorCode ierr;
1972   DM             dm;
1973   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
1974   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
1975   TSIJacobian    ijac;
1976   TSRHSJacobian  rhsjac;
1977 
1978   PetscFunctionBegin;
1979   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1980   if (ts->setupcalled) PetscFunctionReturn(0);
1981 
1982   ts->total_steps = 0;
1983   if (!((PetscObject)ts)->type_name) {
1984     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1985   }
1986 
1987   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1988 
1989 
1990   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
1991 
1992   if (ts->rhsjacobian.reuse) {
1993     Mat Amat,Pmat;
1994     SNES snes;
1995     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1996     ierr = SNESGetJacobian(snes,&Amat,&Pmat,NULL,NULL);CHKERRQ(ierr);
1997     /* Matching matrices implies that an IJacobian is NOT set, because if it had been set, the IJacobian's matrix would
1998      * have displaced the RHS matrix */
1999     if (Amat == ts->Arhs) {
2000       ierr = MatDuplicate(ts->Arhs,MAT_DO_NOT_COPY_VALUES,&Amat);CHKERRQ(ierr);
2001       ierr = SNESSetJacobian(snes,Amat,NULL,NULL,NULL);CHKERRQ(ierr);
2002       ierr = MatDestroy(&Amat);CHKERRQ(ierr);
2003     }
2004     if (Pmat == ts->Brhs) {
2005       ierr = MatDuplicate(ts->Brhs,MAT_DO_NOT_COPY_VALUES,&Pmat);CHKERRQ(ierr);
2006       ierr = SNESSetJacobian(snes,NULL,Pmat,NULL,NULL);CHKERRQ(ierr);
2007       ierr = MatDestroy(&Pmat);CHKERRQ(ierr);
2008     }
2009   }
2010   if (ts->ops->setup) {
2011     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
2012   }
2013 
2014   /* in the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
2015    to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
2016    */
2017   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2018   ierr = DMSNESGetFunction(dm,&func,NULL);CHKERRQ(ierr);
2019   if (!func) {
2020     ierr =DMSNESSetFunction(dm,SNESTSFormFunction,ts);CHKERRQ(ierr);
2021   }
2022   /* if the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
2023      Otherwise, the SNES will use coloring internally to form the Jacobian.
2024    */
2025   ierr = DMSNESGetJacobian(dm,&jac,NULL);CHKERRQ(ierr);
2026   ierr = DMTSGetIJacobian(dm,&ijac,NULL);CHKERRQ(ierr);
2027   ierr = DMTSGetRHSJacobian(dm,&rhsjac,NULL);CHKERRQ(ierr);
2028   if (!jac && (ijac || rhsjac)) {
2029     ierr = DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);CHKERRQ(ierr);
2030   }
2031   ts->setupcalled = PETSC_TRUE;
2032   PetscFunctionReturn(0);
2033 }
2034 
2035 #undef __FUNCT__
2036 #define __FUNCT__ "TSAdjointSetUp"
2037 /*@
2038    TSAdjointSetUp - Sets up the internal data structures for the later use
2039    of an adjoint solver
2040 
2041    Collective on TS
2042 
2043    Input Parameter:
2044 .  ts - the TS context obtained from TSCreate()
2045 
2046    Level: advanced
2047 
2048 .keywords: TS, timestep, setup
2049 
2050 .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients()
2051 @*/
2052 PetscErrorCode  TSAdjointSetUp(TS ts)
2053 {
2054   PetscErrorCode ierr;
2055 
2056   PetscFunctionBegin;
2057   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2058   if (ts->adjointsetupcalled) PetscFunctionReturn(0);
2059   if (!ts->vecs_sensi) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first");
2060 
2061   if (ts->vec_costintegral) { /* if there is integral in the cost function*/
2062     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&ts->vecs_drdy);CHKERRQ(ierr);
2063     if (ts->vecs_sensip){
2064       ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr);
2065     }
2066   }
2067 
2068   if (ts->ops->adjointsetup) {
2069     ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr);
2070   }
2071   ts->adjointsetupcalled = PETSC_TRUE;
2072   PetscFunctionReturn(0);
2073 }
2074 
2075 #undef __FUNCT__
2076 #define __FUNCT__ "TSReset"
2077 /*@
2078    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
2079 
2080    Collective on TS
2081 
2082    Input Parameter:
2083 .  ts - the TS context obtained from TSCreate()
2084 
2085    Level: beginner
2086 
2087 .keywords: TS, timestep, reset
2088 
2089 .seealso: TSCreate(), TSSetup(), TSDestroy()
2090 @*/
2091 PetscErrorCode  TSReset(TS ts)
2092 {
2093   PetscErrorCode ierr;
2094 
2095   PetscFunctionBegin;
2096   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2097 
2098   if (ts->ops->reset) {
2099     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
2100   }
2101   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
2102   if (ts->adapt) {ierr = TSAdaptReset(ts->adapt);CHKERRQ(ierr);}
2103 
2104   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
2105   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
2106   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
2107   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
2108   ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
2109   ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
2110   ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);
2111 
2112  if (ts->vec_costintegral) {
2113     ierr = VecDestroyVecs(ts->numcost,&ts->vecs_drdy);CHKERRQ(ierr);
2114     if (ts->vecs_drdp){
2115       ierr = VecDestroyVecs(ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr);
2116     }
2117   }
2118   ts->vecs_sensi  = NULL;
2119   ts->vecs_sensip = NULL;
2120   ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
2121   ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr);
2122   ierr = VecDestroy(&ts->vec_costintegrand);CHKERRQ(ierr);
2123   ts->setupcalled = PETSC_FALSE;
2124   PetscFunctionReturn(0);
2125 }
2126 
2127 #undef __FUNCT__
2128 #define __FUNCT__ "TSDestroy"
2129 /*@
2130    TSDestroy - Destroys the timestepper context that was created
2131    with TSCreate().
2132 
2133    Collective on TS
2134 
2135    Input Parameter:
2136 .  ts - the TS context obtained from TSCreate()
2137 
2138    Level: beginner
2139 
2140 .keywords: TS, timestepper, destroy
2141 
2142 .seealso: TSCreate(), TSSetUp(), TSSolve()
2143 @*/
2144 PetscErrorCode  TSDestroy(TS *ts)
2145 {
2146   PetscErrorCode ierr;
2147 
2148   PetscFunctionBegin;
2149   if (!*ts) PetscFunctionReturn(0);
2150   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
2151   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
2152 
2153   ierr = TSReset((*ts));CHKERRQ(ierr);
2154 
2155   /* if memory was published with SAWs then destroy it */
2156   ierr = PetscObjectSAWsViewOff((PetscObject)*ts);CHKERRQ(ierr);
2157   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
2158 
2159   ierr = TSTrajectoryDestroy(&(*ts)->trajectory);CHKERRQ(ierr);
2160 
2161   ierr = TSAdaptDestroy(&(*ts)->adapt);CHKERRQ(ierr);
2162   if ((*ts)->event) {
2163     ierr = TSEventMonitorDestroy(&(*ts)->event);CHKERRQ(ierr);
2164   }
2165   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
2166   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
2167   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
2168   ierr = TSAdjointMonitorCancel((*ts));CHKERRQ(ierr);
2169 
2170   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
2171   PetscFunctionReturn(0);
2172 }
2173 
2174 #undef __FUNCT__
2175 #define __FUNCT__ "TSGetSNES"
2176 /*@
2177    TSGetSNES - Returns the SNES (nonlinear solver) associated with
2178    a TS (timestepper) context. Valid only for nonlinear problems.
2179 
2180    Not Collective, but SNES is parallel if TS is parallel
2181 
2182    Input Parameter:
2183 .  ts - the TS context obtained from TSCreate()
2184 
2185    Output Parameter:
2186 .  snes - the nonlinear solver context
2187 
2188    Notes:
2189    The user can then directly manipulate the SNES context to set various
2190    options, etc.  Likewise, the user can then extract and manipulate the
2191    KSP, KSP, and PC contexts as well.
2192 
2193    TSGetSNES() does not work for integrators that do not use SNES; in
2194    this case TSGetSNES() returns NULL in snes.
2195 
2196    Level: beginner
2197 
2198 .keywords: timestep, get, SNES
2199 @*/
2200 PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
2201 {
2202   PetscErrorCode ierr;
2203 
2204   PetscFunctionBegin;
2205   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2206   PetscValidPointer(snes,2);
2207   if (!ts->snes) {
2208     ierr = SNESCreate(PetscObjectComm((PetscObject)ts),&ts->snes);CHKERRQ(ierr);
2209     ierr = SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);CHKERRQ(ierr);
2210     ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->snes);CHKERRQ(ierr);
2211     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
2212     if (ts->dm) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
2213     if (ts->problem_type == TS_LINEAR) {
2214       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
2215     }
2216   }
2217   *snes = ts->snes;
2218   PetscFunctionReturn(0);
2219 }
2220 
2221 #undef __FUNCT__
2222 #define __FUNCT__ "TSSetSNES"
2223 /*@
2224    TSSetSNES - Set the SNES (nonlinear solver) to be used by the timestepping context
2225 
2226    Collective
2227 
2228    Input Parameter:
2229 +  ts - the TS context obtained from TSCreate()
2230 -  snes - the nonlinear solver context
2231 
2232    Notes:
2233    Most users should have the TS created by calling TSGetSNES()
2234 
2235    Level: developer
2236 
2237 .keywords: timestep, set, SNES
2238 @*/
2239 PetscErrorCode TSSetSNES(TS ts,SNES snes)
2240 {
2241   PetscErrorCode ierr;
2242   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
2243 
2244   PetscFunctionBegin;
2245   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2246   PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
2247   ierr = PetscObjectReference((PetscObject)snes);CHKERRQ(ierr);
2248   ierr = SNESDestroy(&ts->snes);CHKERRQ(ierr);
2249 
2250   ts->snes = snes;
2251 
2252   ierr = SNESSetFunction(ts->snes,NULL,SNESTSFormFunction,ts);CHKERRQ(ierr);
2253   ierr = SNESGetJacobian(ts->snes,NULL,NULL,&func,NULL);CHKERRQ(ierr);
2254   if (func == SNESTSFormJacobian) {
2255     ierr = SNESSetJacobian(ts->snes,NULL,NULL,SNESTSFormJacobian,ts);CHKERRQ(ierr);
2256   }
2257   PetscFunctionReturn(0);
2258 }
2259 
2260 #undef __FUNCT__
2261 #define __FUNCT__ "TSGetKSP"
2262 /*@
2263    TSGetKSP - Returns the KSP (linear solver) associated with
2264    a TS (timestepper) context.
2265 
2266    Not Collective, but KSP is parallel if TS is parallel
2267 
2268    Input Parameter:
2269 .  ts - the TS context obtained from TSCreate()
2270 
2271    Output Parameter:
2272 .  ksp - the nonlinear solver context
2273 
2274    Notes:
2275    The user can then directly manipulate the KSP context to set various
2276    options, etc.  Likewise, the user can then extract and manipulate the
2277    KSP and PC contexts as well.
2278 
2279    TSGetKSP() does not work for integrators that do not use KSP;
2280    in this case TSGetKSP() returns NULL in ksp.
2281 
2282    Level: beginner
2283 
2284 .keywords: timestep, get, KSP
2285 @*/
2286 PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
2287 {
2288   PetscErrorCode ierr;
2289   SNES           snes;
2290 
2291   PetscFunctionBegin;
2292   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2293   PetscValidPointer(ksp,2);
2294   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
2295   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
2296   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2297   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
2298   PetscFunctionReturn(0);
2299 }
2300 
2301 /* ----------- Routines to set solver parameters ---------- */
2302 
2303 #undef __FUNCT__
2304 #define __FUNCT__ "TSGetDuration"
2305 /*@
2306    TSGetDuration - Gets the maximum number of timesteps to use and
2307    maximum time for iteration.
2308 
2309    Not Collective
2310 
2311    Input Parameters:
2312 +  ts       - the TS context obtained from TSCreate()
2313 .  maxsteps - maximum number of iterations to use, or NULL
2314 -  maxtime  - final time to iterate to, or NULL
2315 
2316    Level: intermediate
2317 
2318 .keywords: TS, timestep, get, maximum, iterations, time
2319 @*/
2320 PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
2321 {
2322   PetscFunctionBegin;
2323   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2324   if (maxsteps) {
2325     PetscValidIntPointer(maxsteps,2);
2326     *maxsteps = ts->max_steps;
2327   }
2328   if (maxtime) {
2329     PetscValidScalarPointer(maxtime,3);
2330     *maxtime = ts->max_time;
2331   }
2332   PetscFunctionReturn(0);
2333 }
2334 
2335 #undef __FUNCT__
2336 #define __FUNCT__ "TSSetDuration"
2337 /*@
2338    TSSetDuration - Sets the maximum number of timesteps to use and
2339    maximum time for iteration.
2340 
2341    Logically Collective on TS
2342 
2343    Input Parameters:
2344 +  ts - the TS context obtained from TSCreate()
2345 .  maxsteps - maximum number of iterations to use
2346 -  maxtime - final time to iterate to
2347 
2348    Options Database Keys:
2349 .  -ts_max_steps <maxsteps> - Sets maxsteps
2350 .  -ts_final_time <maxtime> - Sets maxtime
2351 
2352    Notes:
2353    The default maximum number of iterations is 5000. Default time is 5.0
2354 
2355    Level: intermediate
2356 
2357 .keywords: TS, timestep, set, maximum, iterations
2358 
2359 .seealso: TSSetExactFinalTime()
2360 @*/
2361 PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
2362 {
2363   PetscFunctionBegin;
2364   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2365   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
2366   PetscValidLogicalCollectiveReal(ts,maxtime,2);
2367   if (maxsteps >= 0) ts->max_steps = maxsteps;
2368   if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime;
2369   PetscFunctionReturn(0);
2370 }
2371 
2372 #undef __FUNCT__
2373 #define __FUNCT__ "TSSetSolution"
2374 /*@
2375    TSSetSolution - Sets the initial solution vector
2376    for use by the TS routines.
2377 
2378    Logically Collective on TS and Vec
2379 
2380    Input Parameters:
2381 +  ts - the TS context obtained from TSCreate()
2382 -  u - the solution vector
2383 
2384    Level: beginner
2385 
2386 .keywords: TS, timestep, set, solution, initial conditions
2387 @*/
2388 PetscErrorCode  TSSetSolution(TS ts,Vec u)
2389 {
2390   PetscErrorCode ierr;
2391   DM             dm;
2392 
2393   PetscFunctionBegin;
2394   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2395   PetscValidHeaderSpecific(u,VEC_CLASSID,2);
2396   ierr = PetscObjectReference((PetscObject)u);CHKERRQ(ierr);
2397   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
2398 
2399   ts->vec_sol = u;
2400 
2401   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2402   ierr = DMShellSetGlobalVector(dm,u);CHKERRQ(ierr);
2403   PetscFunctionReturn(0);
2404 }
2405 
2406 #undef __FUNCT__
2407 #define __FUNCT__ "TSAdjointSetSteps"
2408 /*@
2409    TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
2410 
2411    Logically Collective on TS
2412 
2413    Input Parameters:
2414 +  ts - the TS context obtained from TSCreate()
2415 .  steps - number of steps to use
2416 
2417    Level: intermediate
2418 
2419    Notes: Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this
2420           so as to integrate back to less than the original timestep
2421 
2422 .keywords: TS, timestep, set, maximum, iterations
2423 
2424 .seealso: TSSetExactFinalTime()
2425 @*/
2426 PetscErrorCode  TSAdjointSetSteps(TS ts,PetscInt steps)
2427 {
2428   PetscFunctionBegin;
2429   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2430   PetscValidLogicalCollectiveInt(ts,steps,2);
2431   if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps");
2432   if (steps > ts->total_steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps");
2433   ts->adjoint_max_steps = steps;
2434   PetscFunctionReturn(0);
2435 }
2436 
2437 #undef __FUNCT__
2438 #define __FUNCT__ "TSSetCostGradients"
2439 /*@
2440    TSSetCostGradients - Sets the initial value of the gradients of the cost function w.r.t. initial conditions and w.r.t. the problem parameters
2441       for use by the TSAdjoint routines.
2442 
2443    Logically Collective on TS and Vec
2444 
2445    Input Parameters:
2446 +  ts - the TS context obtained from TSCreate()
2447 .  lambda - gradients with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector
2448 -  mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
2449 
2450    Level: beginner
2451 
2452    Notes: the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime  mu_i = df/dp|finaltime
2453 
2454 .keywords: TS, timestep, set, sensitivity, initial conditions
2455 @*/
2456 PetscErrorCode  TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu)
2457 {
2458   PetscFunctionBegin;
2459   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2460   PetscValidPointer(lambda,2);
2461   ts->vecs_sensi  = lambda;
2462   ts->vecs_sensip = mu;
2463   if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand");
2464   ts->numcost  = numcost;
2465   PetscFunctionReturn(0);
2466 }
2467 
2468 #undef __FUNCT__
2469 #define __FUNCT__ "TSAdjointSetRHSJacobian"
2470 /*@C
2471   TSAdjointSetRHSJacobian - Sets the function that computes the Jacobian of G w.r.t. the parameters p where y_t = G(y,p,t), as well as the location to store the matrix.
2472 
2473   Logically Collective on TS
2474 
2475   Input Parameters:
2476 + ts   - The TS context obtained from TSCreate()
2477 - func - The function
2478 
2479   Calling sequence of func:
2480 $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
2481 +   t - current timestep
2482 .   y - input vector (current ODE solution)
2483 .   A - output matrix
2484 -   ctx - [optional] user-defined function context
2485 
2486   Level: intermediate
2487 
2488   Notes: Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p
2489 
2490 .keywords: TS, sensitivity
2491 .seealso:
2492 @*/
2493 PetscErrorCode  TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
2494 {
2495   PetscErrorCode ierr;
2496 
2497   PetscFunctionBegin;
2498   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2499   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
2500 
2501   ts->rhsjacobianp    = func;
2502   ts->rhsjacobianpctx = ctx;
2503   if(Amat) {
2504     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
2505     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
2506     ts->Jacp = Amat;
2507   }
2508   PetscFunctionReturn(0);
2509 }
2510 
2511 #undef __FUNCT__
2512 #define __FUNCT__ "TSAdjointComputeRHSJacobian"
2513 /*@C
2514   TSAdjointComputeRHSJacobian - Runs the user-defined Jacobian function.
2515 
2516   Collective on TS
2517 
2518   Input Parameters:
2519 . ts   - The TS context obtained from TSCreate()
2520 
2521   Level: developer
2522 
2523 .keywords: TS, sensitivity
2524 .seealso: TSAdjointSetRHSJacobian()
2525 @*/
2526 PetscErrorCode  TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat)
2527 {
2528   PetscErrorCode ierr;
2529 
2530   PetscFunctionBegin;
2531   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2532   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
2533   PetscValidPointer(Amat,4);
2534 
2535   PetscStackPush("TS user JacobianP function for sensitivity analysis");
2536   ierr = (*ts->rhsjacobianp)(ts,t,X,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr);
2537   PetscStackPop;
2538   PetscFunctionReturn(0);
2539 }
2540 
2541 #undef __FUNCT__
2542 #define __FUNCT__ "TSSetCostIntegrand"
2543 /*@C
2544     TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
2545 
2546     Logically Collective on TS
2547 
2548     Input Parameters:
2549 +   ts - the TS context obtained from TSCreate()
2550 .   numcost - number of gradients to be computed, this is the number of cost functions
2551 .   rf - routine for evaluating the integrand function
2552 .   drdyf - function that computes the gradients of the r's with respect to y,NULL if not a function y
2553 .   drdpf - function that computes the gradients of the r's with respect to p, NULL if not a function of p
2554 -   ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
2555 
2556     Calling sequence of rf:
2557 $     rf(TS ts,PetscReal t,Vec y,Vec f[],void *ctx);
2558 
2559 +   t - current timestep
2560 .   y - input vector
2561 .   f - function result; one vector entry for each cost function
2562 -   ctx - [optional] user-defined function context
2563 
2564    Calling sequence of drdyf:
2565 $    PetscErroCode drdyf(TS ts,PetscReal t,Vec y,Vec *drdy,void *ctx);
2566 
2567    Calling sequence of drdpf:
2568 $    PetscErroCode drdpf(TS ts,PetscReal t,Vec y,Vec *drdp,void *ctx);
2569 
2570     Level: intermediate
2571 
2572     Notes: For optimization there is generally a single cost function, numcost = 1. For sensitivities there may be multiple cost functions
2573 
2574 .keywords: TS, sensitivity analysis, timestep, set, quadrature, function
2575 
2576 .seealso: TSAdjointSetRHSJacobian(),TSGetCostGradients(), TSSetCostGradients()
2577 @*/
2578 PetscErrorCode  TSSetCostIntegrand(TS ts,PetscInt numcost, PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
2579                                                                   PetscErrorCode (*drdyf)(TS,PetscReal,Vec,Vec*,void*),
2580                                                                   PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),void *ctx)
2581 {
2582   PetscErrorCode ierr;
2583 
2584   PetscFunctionBegin;
2585   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2586   if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostGradients()");
2587   if (!ts->numcost) ts->numcost=numcost;
2588 
2589   ierr                 = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr);
2590   ierr                 = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr);
2591   ts->costintegrand    = rf;
2592   ts->costintegrandctx = ctx;
2593   ts->drdyfunction     = drdyf;
2594   ts->drdpfunction     = drdpf;
2595   PetscFunctionReturn(0);
2596 }
2597 
2598 #undef __FUNCT__
2599 #define __FUNCT__ "TSGetCostIntegral"
2600 /*@
2601    TSGetCostIntegral - Returns the values of the integral term in the cost functions.
2602    It is valid to call the routine after a backward run.
2603 
2604    Not Collective
2605 
2606    Input Parameter:
2607 .  ts - the TS context obtained from TSCreate()
2608 
2609    Output Parameter:
2610 .  v - the vector containing the integrals for each cost function
2611 
2612    Level: intermediate
2613 
2614 .seealso: TSSetCostIntegrand()
2615 
2616 .keywords: TS, sensitivity analysis
2617 @*/
2618 PetscErrorCode  TSGetCostIntegral(TS ts,Vec *v)
2619 {
2620   PetscFunctionBegin;
2621   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2622   PetscValidPointer(v,2);
2623   *v = ts->vec_costintegral;
2624   PetscFunctionReturn(0);
2625 }
2626 
2627 #undef __FUNCT__
2628 #define __FUNCT__ "TSAdjointComputeCostIntegrand"
2629 /*@
2630    TSAdjointComputeCostIntegrand - Evaluates the integral function in the cost functions.
2631 
2632    Input Parameters:
2633 +  ts - the TS context
2634 .  t - current time
2635 -  y - state vector, i.e. current solution
2636 
2637    Output Parameter:
2638 .  q - vector of size numcost to hold the outputs
2639 
2640    Note:
2641    Most users should not need to explicitly call this routine, as it
2642    is used internally within the sensitivity analysis context.
2643 
2644    Level: developer
2645 
2646 .keywords: TS, compute
2647 
2648 .seealso: TSSetCostIntegrand()
2649 @*/
2650 PetscErrorCode TSAdjointComputeCostIntegrand(TS ts,PetscReal t,Vec y,Vec q)
2651 {
2652   PetscErrorCode ierr;
2653 
2654   PetscFunctionBegin;
2655   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2656   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
2657   PetscValidHeaderSpecific(q,VEC_CLASSID,4);
2658 
2659   ierr = PetscLogEventBegin(TS_FunctionEval,ts,y,q,0);CHKERRQ(ierr);
2660   if (ts->costintegrand) {
2661     PetscStackPush("TS user integrand in the cost function");
2662     ierr = (*ts->costintegrand)(ts,t,y,q,ts->costintegrandctx);CHKERRQ(ierr);
2663     PetscStackPop;
2664   } else {
2665     ierr = VecZeroEntries(q);CHKERRQ(ierr);
2666   }
2667 
2668   ierr = PetscLogEventEnd(TS_FunctionEval,ts,y,q,0);CHKERRQ(ierr);
2669   PetscFunctionReturn(0);
2670 }
2671 
2672 #undef __FUNCT__
2673 #define __FUNCT__ "TSAdjointComputeDRDYFunction"
2674 /*@
2675   TSAdjointComputeDRDYFunction - Runs the user-defined DRDY function.
2676 
2677   Collective on TS
2678 
2679   Input Parameters:
2680 . ts   - The TS context obtained from TSCreate()
2681 
2682   Notes:
2683   TSAdjointComputeDRDYFunction() is typically used for sensitivity implementation,
2684   so most users would not generally call this routine themselves.
2685 
2686   Level: developer
2687 
2688 .keywords: TS, sensitivity
2689 .seealso: TSAdjointComputeDRDYFunction()
2690 @*/
2691 PetscErrorCode  TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec y,Vec *drdy)
2692 {
2693   PetscErrorCode ierr;
2694 
2695   PetscFunctionBegin;
2696   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2697   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
2698 
2699   PetscStackPush("TS user DRDY function for sensitivity analysis");
2700   ierr = (*ts->drdyfunction)(ts,t,y,drdy,ts->costintegrandctx); CHKERRQ(ierr);
2701   PetscStackPop;
2702   PetscFunctionReturn(0);
2703 }
2704 
2705 #undef __FUNCT__
2706 #define __FUNCT__ "TSAdjointComputeDRDPFunction"
2707 /*@
2708   TSAdjointComputeDRDPFunction - Runs the user-defined DRDP function.
2709 
2710   Collective on TS
2711 
2712   Input Parameters:
2713 . ts   - The TS context obtained from TSCreate()
2714 
2715   Notes:
2716   TSDRDPFunction() is typically used for sensitivity implementation,
2717   so most users would not generally call this routine themselves.
2718 
2719   Level: developer
2720 
2721 .keywords: TS, sensitivity
2722 .seealso: TSAdjointSetDRDPFunction()
2723 @*/
2724 PetscErrorCode  TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec y,Vec *drdp)
2725 {
2726   PetscErrorCode ierr;
2727 
2728   PetscFunctionBegin;
2729   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2730   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
2731 
2732   PetscStackPush("TS user DRDP function for sensitivity analysis");
2733   ierr = (*ts->drdpfunction)(ts,t,y,drdp,ts->costintegrandctx); CHKERRQ(ierr);
2734   PetscStackPop;
2735   PetscFunctionReturn(0);
2736 }
2737 
2738 #undef __FUNCT__
2739 #define __FUNCT__ "TSSetPreStep"
2740 /*@C
2741   TSSetPreStep - Sets the general-purpose function
2742   called once at the beginning of each time step.
2743 
2744   Logically Collective on TS
2745 
2746   Input Parameters:
2747 + ts   - The TS context obtained from TSCreate()
2748 - func - The function
2749 
2750   Calling sequence of func:
2751 . func (TS ts);
2752 
2753   Level: intermediate
2754 
2755   Note:
2756   If a step is rejected, TSStep() will call this routine again before each attempt.
2757   The last completed time step number can be queried using TSGetTimeStepNumber(), the
2758   size of the step being attempted can be obtained using TSGetTimeStep().
2759 
2760 .keywords: TS, timestep
2761 .seealso: TSSetPreStage(), TSSetPostStage(), TSSetPostStep(), TSStep()
2762 @*/
2763 PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
2764 {
2765   PetscFunctionBegin;
2766   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2767   ts->prestep = func;
2768   PetscFunctionReturn(0);
2769 }
2770 
2771 #undef __FUNCT__
2772 #define __FUNCT__ "TSPreStep"
2773 /*@
2774   TSPreStep - Runs the user-defined pre-step function.
2775 
2776   Collective on TS
2777 
2778   Input Parameters:
2779 . ts   - The TS context obtained from TSCreate()
2780 
2781   Notes:
2782   TSPreStep() is typically used within time stepping implementations,
2783   so most users would not generally call this routine themselves.
2784 
2785   Level: developer
2786 
2787 .keywords: TS, timestep
2788 .seealso: TSSetPreStep(), TSPreStage(), TSPostStage(), TSPostStep()
2789 @*/
2790 PetscErrorCode  TSPreStep(TS ts)
2791 {
2792   PetscErrorCode ierr;
2793 
2794   PetscFunctionBegin;
2795   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2796   if (ts->prestep) {
2797     PetscStackCallStandard((*ts->prestep),(ts));
2798   }
2799   PetscFunctionReturn(0);
2800 }
2801 
2802 #undef __FUNCT__
2803 #define __FUNCT__ "TSSetPreStage"
2804 /*@C
2805   TSSetPreStage - Sets the general-purpose function
2806   called once at the beginning of each stage.
2807 
2808   Logically Collective on TS
2809 
2810   Input Parameters:
2811 + ts   - The TS context obtained from TSCreate()
2812 - func - The function
2813 
2814   Calling sequence of func:
2815 . PetscErrorCode func(TS ts, PetscReal stagetime);
2816 
2817   Level: intermediate
2818 
2819   Note:
2820   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
2821   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
2822   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
2823 
2824 .keywords: TS, timestep
2825 .seealso: TSSetPostStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
2826 @*/
2827 PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
2828 {
2829   PetscFunctionBegin;
2830   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2831   ts->prestage = func;
2832   PetscFunctionReturn(0);
2833 }
2834 
2835 #undef __FUNCT__
2836 #define __FUNCT__ "TSSetPostStage"
2837 /*@C
2838   TSSetPostStage - Sets the general-purpose function
2839   called once at the end of each stage.
2840 
2841   Logically Collective on TS
2842 
2843   Input Parameters:
2844 + ts   - The TS context obtained from TSCreate()
2845 - func - The function
2846 
2847   Calling sequence of func:
2848 . PetscErrorCode func(TS ts, PetscReal stagetime, PetscInt stageindex, Vec* Y);
2849 
2850   Level: intermediate
2851 
2852   Note:
2853   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
2854   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
2855   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
2856 
2857 .keywords: TS, timestep
2858 .seealso: TSSetPreStage(), TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
2859 @*/
2860 PetscErrorCode  TSSetPostStage(TS ts, PetscErrorCode (*func)(TS,PetscReal,PetscInt,Vec*))
2861 {
2862   PetscFunctionBegin;
2863   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2864   ts->poststage = func;
2865   PetscFunctionReturn(0);
2866 }
2867 
2868 #undef __FUNCT__
2869 #define __FUNCT__ "TSPreStage"
2870 /*@
2871   TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()
2872 
2873   Collective on TS
2874 
2875   Input Parameters:
2876 . ts          - The TS context obtained from TSCreate()
2877   stagetime   - The absolute time of the current stage
2878 
2879   Notes:
2880   TSPreStage() is typically used within time stepping implementations,
2881   most users would not generally call this routine themselves.
2882 
2883   Level: developer
2884 
2885 .keywords: TS, timestep
2886 .seealso: TSPostStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
2887 @*/
2888 PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
2889 {
2890   PetscErrorCode ierr;
2891 
2892   PetscFunctionBegin;
2893   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2894   if (ts->prestage) {
2895     PetscStackCallStandard((*ts->prestage),(ts,stagetime));
2896   }
2897   PetscFunctionReturn(0);
2898 }
2899 
2900 #undef __FUNCT__
2901 #define __FUNCT__ "TSPostStage"
2902 /*@
2903   TSPostStage - Runs the user-defined post-stage function set using TSSetPostStage()
2904 
2905   Collective on TS
2906 
2907   Input Parameters:
2908 . ts          - The TS context obtained from TSCreate()
2909   stagetime   - The absolute time of the current stage
2910   stageindex  - Stage number
2911   Y           - Array of vectors (of size = total number
2912                 of stages) with the stage solutions
2913 
2914   Notes:
2915   TSPostStage() is typically used within time stepping implementations,
2916   most users would not generally call this routine themselves.
2917 
2918   Level: developer
2919 
2920 .keywords: TS, timestep
2921 .seealso: TSPreStage(), TSSetPreStep(), TSPreStep(), TSPostStep()
2922 @*/
2923 PetscErrorCode  TSPostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
2924 {
2925   PetscErrorCode ierr;
2926 
2927   PetscFunctionBegin;
2928   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2929   if (ts->poststage) {
2930     PetscStackCallStandard((*ts->poststage),(ts,stagetime,stageindex,Y));
2931   }
2932   PetscFunctionReturn(0);
2933 }
2934 
2935 #undef __FUNCT__
2936 #define __FUNCT__ "TSSetPostStep"
2937 /*@C
2938   TSSetPostStep - Sets the general-purpose function
2939   called once at the end of each time step.
2940 
2941   Logically Collective on TS
2942 
2943   Input Parameters:
2944 + ts   - The TS context obtained from TSCreate()
2945 - func - The function
2946 
2947   Calling sequence of func:
2948 $ func (TS ts);
2949 
2950   Level: intermediate
2951 
2952 .keywords: TS, timestep
2953 .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime()
2954 @*/
2955 PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
2956 {
2957   PetscFunctionBegin;
2958   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2959   ts->poststep = func;
2960   PetscFunctionReturn(0);
2961 }
2962 
2963 #undef __FUNCT__
2964 #define __FUNCT__ "TSPostStep"
2965 /*@
2966   TSPostStep - Runs the user-defined post-step function.
2967 
2968   Collective on TS
2969 
2970   Input Parameters:
2971 . ts   - The TS context obtained from TSCreate()
2972 
2973   Notes:
2974   TSPostStep() is typically used within time stepping implementations,
2975   so most users would not generally call this routine themselves.
2976 
2977   Level: developer
2978 
2979 .keywords: TS, timestep
2980 @*/
2981 PetscErrorCode  TSPostStep(TS ts)
2982 {
2983   PetscErrorCode ierr;
2984 
2985   PetscFunctionBegin;
2986   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2987   if (ts->poststep) {
2988     PetscStackCallStandard((*ts->poststep),(ts));
2989   }
2990   PetscFunctionReturn(0);
2991 }
2992 
2993 /* ------------ Routines to set performance monitoring options ----------- */
2994 
2995 #undef __FUNCT__
2996 #define __FUNCT__ "TSMonitorSet"
2997 /*@C
2998    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
2999    timestep to display the iteration's  progress.
3000 
3001    Logically Collective on TS
3002 
3003    Input Parameters:
3004 +  ts - the TS context obtained from TSCreate()
3005 .  monitor - monitoring routine
3006 .  mctx - [optional] user-defined context for private data for the
3007              monitor routine (use NULL if no context is desired)
3008 -  monitordestroy - [optional] routine that frees monitor context
3009           (may be NULL)
3010 
3011    Calling sequence of monitor:
3012 $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx)
3013 
3014 +    ts - the TS context
3015 .    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
3016                                been interpolated to)
3017 .    time - current time
3018 .    u - current iterate
3019 -    mctx - [optional] monitoring context
3020 
3021    Notes:
3022    This routine adds an additional monitor to the list of monitors that
3023    already has been loaded.
3024 
3025    Fortran notes: Only a single monitor function can be set for each TS object
3026 
3027    Level: intermediate
3028 
3029 .keywords: TS, timestep, set, monitor
3030 
3031 .seealso: TSMonitorDefault(), TSMonitorCancel()
3032 @*/
3033 PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
3034 {
3035   PetscFunctionBegin;
3036   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3037   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3038   ts->monitor[ts->numbermonitors]          = monitor;
3039   ts->monitordestroy[ts->numbermonitors]   = mdestroy;
3040   ts->monitorcontext[ts->numbermonitors++] = (void*)mctx;
3041   PetscFunctionReturn(0);
3042 }
3043 
3044 #undef __FUNCT__
3045 #define __FUNCT__ "TSMonitorCancel"
3046 /*@C
3047    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
3048 
3049    Logically Collective on TS
3050 
3051    Input Parameters:
3052 .  ts - the TS context obtained from TSCreate()
3053 
3054    Notes:
3055    There is no way to remove a single, specific monitor.
3056 
3057    Level: intermediate
3058 
3059 .keywords: TS, timestep, set, monitor
3060 
3061 .seealso: TSMonitorDefault(), TSMonitorSet()
3062 @*/
3063 PetscErrorCode  TSMonitorCancel(TS ts)
3064 {
3065   PetscErrorCode ierr;
3066   PetscInt       i;
3067 
3068   PetscFunctionBegin;
3069   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3070   for (i=0; i<ts->numbermonitors; i++) {
3071     if (ts->monitordestroy[i]) {
3072       ierr = (*ts->monitordestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
3073     }
3074   }
3075   ts->numbermonitors = 0;
3076   PetscFunctionReturn(0);
3077 }
3078 
3079 #undef __FUNCT__
3080 #define __FUNCT__ "TSMonitorDefault"
3081 /*@
3082    TSMonitorDefault - Sets the Default monitor
3083 
3084    Level: intermediate
3085 
3086 .keywords: TS, set, monitor
3087 
3088 .seealso: TSMonitorDefault(), TSMonitorSet()
3089 @*/
3090 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
3091 {
3092   PetscErrorCode ierr;
3093   PetscViewer    viewer =  (PetscViewer) dummy;
3094 
3095   PetscFunctionBegin;
3096   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
3097   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3098   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr);
3099   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3100   PetscFunctionReturn(0);
3101 }
3102 
3103 #undef __FUNCT__
3104 #define __FUNCT__ "TSAdjointMonitorSet"
3105 /*@C
3106    TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
3107    timestep to display the iteration's  progress.
3108 
3109    Logically Collective on TS
3110 
3111    Input Parameters:
3112 +  ts - the TS context obtained from TSCreate()
3113 .  adjointmonitor - monitoring routine
3114 .  adjointmctx - [optional] user-defined context for private data for the
3115              monitor routine (use NULL if no context is desired)
3116 -  adjointmonitordestroy - [optional] routine that frees monitor context
3117           (may be NULL)
3118 
3119    Calling sequence of monitor:
3120 $    int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
3121 
3122 +    ts - the TS context
3123 .    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
3124                                been interpolated to)
3125 .    time - current time
3126 .    u - current iterate
3127 .    numcost - number of cost functionos
3128 .    lambda - sensitivities to initial conditions
3129 .    mu - sensitivities to parameters
3130 -    adjointmctx - [optional] adjoint monitoring context
3131 
3132    Notes:
3133    This routine adds an additional monitor to the list of monitors that
3134    already has been loaded.
3135 
3136    Fortran notes: Only a single monitor function can be set for each TS object
3137 
3138    Level: intermediate
3139 
3140 .keywords: TS, timestep, set, adjoint, monitor
3141 
3142 .seealso: TSAdjointMonitorCancel()
3143 @*/
3144 PetscErrorCode  TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**))
3145 {
3146   PetscFunctionBegin;
3147   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3148   if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set");
3149   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
3150   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
3151   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx;
3152   PetscFunctionReturn(0);
3153 }
3154 
3155 #undef __FUNCT__
3156 #define __FUNCT__ "TSAdjointMonitorCancel"
3157 /*@C
3158    TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
3159 
3160    Logically Collective on TS
3161 
3162    Input Parameters:
3163 .  ts - the TS context obtained from TSCreate()
3164 
3165    Notes:
3166    There is no way to remove a single, specific monitor.
3167 
3168    Level: intermediate
3169 
3170 .keywords: TS, timestep, set, adjoint, monitor
3171 
3172 .seealso: TSAdjointMonitorSet()
3173 @*/
3174 PetscErrorCode  TSAdjointMonitorCancel(TS ts)
3175 {
3176   PetscErrorCode ierr;
3177   PetscInt       i;
3178 
3179   PetscFunctionBegin;
3180   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3181   for (i=0; i<ts->numberadjointmonitors; i++) {
3182     if (ts->adjointmonitordestroy[i]) {
3183       ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
3184     }
3185   }
3186   ts->numberadjointmonitors = 0;
3187   PetscFunctionReturn(0);
3188 }
3189 
3190 #undef __FUNCT__
3191 #define __FUNCT__ "TSAdjointMonitorDefault"
3192 /*@
3193    TSAdjointMonitorDefault - Sets the Default monitor
3194 
3195    Level: intermediate
3196 
3197 .keywords: TS, set, monitor
3198 
3199 .seealso: TSAdjointMonitorSet()
3200 @*/
3201 PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
3202 {
3203   PetscErrorCode ierr;
3204   PetscViewer    viewer =  (PetscViewer) dummy;
3205 
3206   PetscFunctionBegin;
3207   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
3208   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3209   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr);
3210   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
3211   PetscFunctionReturn(0);
3212 }
3213 
3214 #undef __FUNCT__
3215 #define __FUNCT__ "TSSetRetainStages"
3216 /*@
3217    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
3218 
3219    Logically Collective on TS
3220 
3221    Input Argument:
3222 .  ts - time stepping context
3223 
3224    Output Argument:
3225 .  flg - PETSC_TRUE or PETSC_FALSE
3226 
3227    Level: intermediate
3228 
3229 .keywords: TS, set
3230 
3231 .seealso: TSInterpolate(), TSSetPostStep()
3232 @*/
3233 PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
3234 {
3235   PetscFunctionBegin;
3236   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3237   ts->retain_stages = flg;
3238   PetscFunctionReturn(0);
3239 }
3240 
3241 #undef __FUNCT__
3242 #define __FUNCT__ "TSInterpolate"
3243 /*@
3244    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
3245 
3246    Collective on TS
3247 
3248    Input Argument:
3249 +  ts - time stepping context
3250 -  t - time to interpolate to
3251 
3252    Output Argument:
3253 .  U - state at given time
3254 
3255    Notes:
3256    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
3257 
3258    Level: intermediate
3259 
3260    Developer Notes:
3261    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
3262 
3263 .keywords: TS, set
3264 
3265 .seealso: TSSetRetainStages(), TSSetPostStep()
3266 @*/
3267 PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec U)
3268 {
3269   PetscErrorCode ierr;
3270 
3271   PetscFunctionBegin;
3272   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3273   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
3274   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,(double)(ts->ptime-ts->time_step_prev),(double)ts->ptime);
3275   if (!ts->ops->interpolate) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
3276   ierr = (*ts->ops->interpolate)(ts,t,U);CHKERRQ(ierr);
3277   PetscFunctionReturn(0);
3278 }
3279 
3280 #undef __FUNCT__
3281 #define __FUNCT__ "TSStep"
3282 /*@
3283    TSStep - Steps one time step
3284 
3285    Collective on TS
3286 
3287    Input Parameter:
3288 .  ts - the TS context obtained from TSCreate()
3289 
3290    Level: developer
3291 
3292    Notes:
3293    The public interface for the ODE/DAE solvers is TSSolve(), you should almost for sure be using that routine and not this routine.
3294 
3295    The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
3296    be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.
3297 
3298    This may over-step the final time provided in TSSetDuration() depending on the time-step used. TSSolve() interpolates to exactly the
3299    time provided in TSSetDuration(). One can use TSInterpolate() to determine an interpolated solution within the final timestep.
3300 
3301 .keywords: TS, timestep, solve
3302 
3303 .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSSetPostStage(), TSInterpolate()
3304 @*/
3305 PetscErrorCode  TSStep(TS ts)
3306 {
3307   DM               dm;
3308   PetscErrorCode   ierr;
3309   static PetscBool cite = PETSC_FALSE;
3310 
3311   PetscFunctionBegin;
3312   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
3313   ierr = PetscCitationsRegister("@techreport{tspaper,\n"
3314                                 "  title       = {{PETSc/TS}: A Modern Scalable {DAE/ODE} Solver Library},\n"
3315                                 "  author      = {Shrirang Abhyankar and Jed Brown and Emil Constantinescu and Debojyoti Ghosh and Barry F. Smith},\n"
3316                                 "  type        = {Preprint},\n"
3317                                 "  number      = {ANL/MCS-P5061-0114},\n"
3318                                 "  institution = {Argonne National Laboratory},\n"
3319                                 "  year        = {2014}\n}\n",&cite);CHKERRQ(ierr);
3320 
3321   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
3322   ierr = TSSetUp(ts);CHKERRQ(ierr);
3323 
3324   ts->reason = TS_CONVERGED_ITERATING;
3325   ts->ptime_prev = ts->ptime;
3326   ierr = DMSetOutputSequenceNumber(dm, ts->steps, ts->ptime);CHKERRQ(ierr);
3327 
3328   if (!ts->ops->step) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSStep not implemented for type '%s'",((PetscObject)ts)->type_name);
3329   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
3330   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
3331   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
3332 
3333   ts->time_step_prev = ts->ptime - ts->ptime_prev;
3334   ierr = DMSetOutputSequenceNumber(dm, ts->steps, ts->ptime);CHKERRQ(ierr);
3335 
3336   if (ts->reason < 0) {
3337     if (ts->errorifstepfailed) {
3338       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) 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]);
3339       else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
3340     }
3341   } else if (!ts->reason) {
3342     if (ts->steps >= ts->max_steps)     ts->reason = TS_CONVERGED_ITS;
3343     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
3344   }
3345   ts->total_steps++;
3346   ts->steprollback = PETSC_FALSE;
3347   PetscFunctionReturn(0);
3348 }
3349 
3350 #undef __FUNCT__
3351 #define __FUNCT__ "TSAdjointStep"
3352 /*@
3353    TSAdjointStep - Steps one time step backward in the adjoint run
3354 
3355    Collective on TS
3356 
3357    Input Parameter:
3358 .  ts - the TS context obtained from TSCreate()
3359 
3360    Level: intermediate
3361 
3362 .keywords: TS, adjoint, step
3363 
3364 .seealso: TSAdjointSetUp(), TSAdjointSolve()
3365 @*/
3366 PetscErrorCode  TSAdjointStep(TS ts)
3367 {
3368   DM               dm;
3369   PetscErrorCode   ierr;
3370 
3371   PetscFunctionBegin;
3372   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
3373   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
3374   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
3375 
3376   ts->reason = TS_CONVERGED_ITERATING;
3377   ts->ptime_prev = ts->ptime;
3378   ierr = DMSetOutputSequenceNumber(dm, ts->steps, ts->ptime);CHKERRQ(ierr);
3379   ierr = VecViewFromOptions(ts->vec_sol,(PetscObject)ts, "-ts_view_solution");CHKERRQ(ierr);
3380 
3381   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
3382   if (!ts->ops->adjointstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed because the adjoint of  %s has not been implemented, try other time stepping methods for adjoint sensitivity analysis",((PetscObject)ts)->type_name);
3383   ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr);
3384   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
3385 
3386   ts->time_step_prev = ts->ptime - ts->ptime_prev;
3387   ierr = DMSetOutputSequenceNumber(dm, ts->steps, ts->ptime);CHKERRQ(ierr);
3388 
3389   if (ts->reason < 0) {
3390     if (ts->errorifstepfailed) {
3391       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) 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]);
3392       else if (ts->reason == TS_DIVERGED_STEP_REJECTED) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_reject or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
3393       else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
3394     }
3395   } else if (!ts->reason) {
3396     if (ts->steps >= ts->adjoint_max_steps)     ts->reason = TS_CONVERGED_ITS;
3397     else if (ts->ptime >= ts->max_time)         ts->reason = TS_CONVERGED_TIME;
3398   }
3399   ts->total_steps--;
3400   PetscFunctionReturn(0);
3401 }
3402 
3403 #undef __FUNCT__
3404 #define __FUNCT__ "TSEvaluateStep"
3405 /*@
3406    TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
3407 
3408    Collective on TS
3409 
3410    Input Arguments:
3411 +  ts - time stepping context
3412 .  order - desired order of accuracy
3413 -  done - whether the step was evaluated at this order (pass NULL to generate an error if not available)
3414 
3415    Output Arguments:
3416 .  U - state at the end of the current step
3417 
3418    Level: advanced
3419 
3420    Notes:
3421    This function cannot be called until all stages have been evaluated.
3422    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.
3423 
3424 .seealso: TSStep(), TSAdapt
3425 @*/
3426 PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec U,PetscBool *done)
3427 {
3428   PetscErrorCode ierr;
3429 
3430   PetscFunctionBegin;
3431   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3432   PetscValidType(ts,1);
3433   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
3434   if (!ts->ops->evaluatestep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
3435   ierr = (*ts->ops->evaluatestep)(ts,order,U,done);CHKERRQ(ierr);
3436   PetscFunctionReturn(0);
3437 }
3438 
3439 
3440 #undef __FUNCT__
3441 #define __FUNCT__ "TSSolve"
3442 /*@
3443    TSSolve - Steps the requested number of timesteps.
3444 
3445    Collective on TS
3446 
3447    Input Parameter:
3448 +  ts - the TS context obtained from TSCreate()
3449 -  u - the solution vector  (can be null if TSSetSolution() was used, otherwise must contain the initial conditions)
3450 
3451    Level: beginner
3452 
3453    Notes:
3454    The final time returned by this function may be different from the time of the internally
3455    held state accessible by TSGetSolution() and TSGetTime() because the method may have
3456    stepped over the final time.
3457 
3458 .keywords: TS, timestep, solve
3459 
3460 .seealso: TSCreate(), TSSetSolution(), TSStep()
3461 @*/
3462 PetscErrorCode TSSolve(TS ts,Vec u)
3463 {
3464   Vec               solution;
3465   PetscErrorCode    ierr;
3466 
3467   PetscFunctionBegin;
3468   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3469   if (u) PetscValidHeaderSpecific(u,VEC_CLASSID,2);
3470   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 */
3471     PetscValidHeaderSpecific(u,VEC_CLASSID,2);
3472     if (!ts->vec_sol || u == ts->vec_sol) {
3473       ierr = VecDuplicate(u,&solution);CHKERRQ(ierr);
3474       ierr = TSSetSolution(ts,solution);CHKERRQ(ierr);
3475       ierr = VecDestroy(&solution);CHKERRQ(ierr); /* grant ownership */
3476     }
3477     ierr = VecCopy(u,ts->vec_sol);CHKERRQ(ierr);
3478   } else if (u) {
3479     ierr = TSSetSolution(ts,u);CHKERRQ(ierr);
3480   }
3481   ierr = TSSetUp(ts);CHKERRQ(ierr);
3482   /* reset time step and iteration counters */
3483   ts->steps             = 0;
3484   ts->ksp_its           = 0;
3485   ts->snes_its          = 0;
3486   ts->num_snes_failures = 0;
3487   ts->reject            = 0;
3488   ts->reason            = TS_CONVERGED_ITERATING;
3489 
3490   ierr = TSViewFromOptions(ts,NULL,"-ts_view_pre");CHKERRQ(ierr);
3491   {
3492     DM dm;
3493     ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
3494     ierr = DMSetOutputSequenceNumber(dm, ts->steps, ts->ptime);CHKERRQ(ierr);
3495   }
3496 
3497   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
3498     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
3499     ierr = VecCopy(ts->vec_sol,u);CHKERRQ(ierr);
3500     ts->solvetime = ts->ptime;
3501   } else {
3502     /* steps the requested number of timesteps. */
3503     if (ts->steps >= ts->max_steps)     ts->reason = TS_CONVERGED_ITS;
3504     else if (ts->ptime >= ts->max_time) ts->reason = TS_CONVERGED_TIME;
3505     ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
3506     if (ts->vec_costintegral) ts->costintegralfwd=PETSC_TRUE;
3507     if(ts->event) {
3508       ierr = TSEventMonitorInitialize(ts);CHKERRQ(ierr);
3509     }
3510     while (!ts->reason) {
3511       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
3512       ierr = TSStep(ts);CHKERRQ(ierr);
3513       if (ts->event) {
3514 	ierr = TSEventMonitor(ts);CHKERRQ(ierr);
3515       }
3516       if(!ts->steprollback) {
3517 	ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
3518 	ierr = TSPostStep(ts);CHKERRQ(ierr);
3519       }
3520     }
3521     if (ts->exact_final_time == TS_EXACTFINALTIME_INTERPOLATE && ts->ptime > ts->max_time) {
3522       ierr = TSInterpolate(ts,ts->max_time,u);CHKERRQ(ierr);
3523       ts->solvetime = ts->max_time;
3524       solution = u;
3525     } else {
3526       if (u) {ierr = VecCopy(ts->vec_sol,u);CHKERRQ(ierr);}
3527       ts->solvetime = ts->ptime;
3528       solution = ts->vec_sol;
3529     }
3530     ierr = TSMonitor(ts,ts->steps,ts->solvetime,solution);CHKERRQ(ierr);
3531     ierr = VecViewFromOptions(solution,(PetscObject) ts,"-ts_view_solution");CHKERRQ(ierr);
3532   }
3533 
3534   ierr = TSViewFromOptions(ts,NULL,"-ts_view");CHKERRQ(ierr);
3535   ierr = VecViewFromOptions(ts->vec_sol,NULL,"-ts_view_solution");CHKERRQ(ierr);
3536   ierr = PetscObjectSAWsBlock((PetscObject)ts);CHKERRQ(ierr);
3537   if (ts->adjoint_solve) {
3538     ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
3539   }
3540   PetscFunctionReturn(0);
3541 }
3542 
3543 #undef __FUNCT__
3544 #define __FUNCT__ "TSAdjointSolve"
3545 /*@
3546    TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
3547 
3548    Collective on TS
3549 
3550    Input Parameter:
3551 .  ts - the TS context obtained from TSCreate()
3552 
3553    Options Database:
3554 . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial conditions
3555 
3556    Level: intermediate
3557 
3558    Notes:
3559    This must be called after a call to TSSolve() that solves the forward problem
3560 
3561    By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time
3562 
3563 .keywords: TS, timestep, solve
3564 
3565 .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep()
3566 @*/
3567 PetscErrorCode TSAdjointSolve(TS ts)
3568 {
3569   PetscErrorCode    ierr;
3570 
3571   PetscFunctionBegin;
3572   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3573   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
3574   /* reset time step and iteration counters */
3575   ts->steps             = 0;
3576   ts->ksp_its           = 0;
3577   ts->snes_its          = 0;
3578   ts->num_snes_failures = 0;
3579   ts->reject            = 0;
3580   ts->reason            = TS_CONVERGED_ITERATING;
3581 
3582   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->total_steps;
3583 
3584   if (ts->steps >= ts->adjoint_max_steps)     ts->reason = TS_CONVERGED_ITS;
3585   while (!ts->reason) {
3586     ierr = TSTrajectoryGet(ts->trajectory,ts,ts->total_steps,&ts->ptime);CHKERRQ(ierr);
3587     ierr = TSAdjointMonitor(ts,ts->total_steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
3588     if (ts->event) {
3589       ierr = TSAdjointEventMonitor(ts);CHKERRQ(ierr);
3590     }
3591     ierr = TSAdjointStep(ts);CHKERRQ(ierr);
3592   }
3593   ierr = TSTrajectoryGet(ts->trajectory,ts,ts->total_steps,&ts->ptime);CHKERRQ(ierr);
3594   ierr = TSAdjointMonitor(ts,ts->total_steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
3595   ts->solvetime = ts->ptime;
3596   ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr);
3597   PetscFunctionReturn(0);
3598 }
3599 
3600 #undef __FUNCT__
3601 #define __FUNCT__ "TSMonitor"
3602 /*@C
3603    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
3604 
3605    Collective on TS
3606 
3607    Input Parameters:
3608 +  ts - time stepping context obtained from TSCreate()
3609 .  step - step number that has just completed
3610 .  ptime - model time of the state
3611 -  u - state at the current model time
3612 
3613    Notes:
3614    TSMonitor() is typically used automatically within the time stepping implementations.
3615    Users would almost never call this routine directly.
3616 
3617    Level: developer
3618 
3619 .keywords: TS, timestep
3620 @*/
3621 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u)
3622 {
3623   PetscErrorCode ierr;
3624   PetscInt       i,n = ts->numbermonitors;
3625 
3626   PetscFunctionBegin;
3627   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3628   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
3629   ierr = VecLockPush(u);CHKERRQ(ierr);
3630   for (i=0; i<n; i++) {
3631     ierr = (*ts->monitor[i])(ts,step,ptime,u,ts->monitorcontext[i]);CHKERRQ(ierr);
3632   }
3633   ierr = VecLockPop(u);CHKERRQ(ierr);
3634   PetscFunctionReturn(0);
3635 }
3636 
3637 #undef __FUNCT__
3638 #define __FUNCT__ "TSAdjointMonitor"
3639 /*@C
3640    TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
3641 
3642    Collective on TS
3643 
3644    Input Parameters:
3645 +  ts - time stepping context obtained from TSCreate()
3646 .  step - step number that has just completed
3647 .  ptime - model time of the state
3648 .  u - state at the current model time
3649 .  numcost - number of cost functions (dimension of lambda  or mu)
3650 .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
3651 -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
3652 
3653    Notes:
3654    TSAdjointMonitor() is typically used automatically within the time stepping implementations.
3655    Users would almost never call this routine directly.
3656 
3657    Level: developer
3658 
3659 .keywords: TS, timestep
3660 @*/
3661 PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu)
3662 {
3663   PetscErrorCode ierr;
3664   PetscInt       i,n = ts->numberadjointmonitors;
3665 
3666   PetscFunctionBegin;
3667   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3668   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
3669   ierr = VecLockPush(u);CHKERRQ(ierr);
3670   for (i=0; i<n; i++) {
3671     ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
3672   }
3673   ierr = VecLockPop(u);CHKERRQ(ierr);
3674   PetscFunctionReturn(0);
3675 }
3676 
3677 /* ------------------------------------------------------------------------*/
3678 #undef __FUNCT__
3679 #define __FUNCT__ "TSMonitorLGCtxCreate"
3680 /*@C
3681    TSMonitorLGCtxCreate - Creates a TSMonitorLGCtx context for use with
3682    TS to monitor the solution process graphically in various ways
3683 
3684    Collective on TS
3685 
3686    Input Parameters:
3687 +  host - the X display to open, or null for the local machine
3688 .  label - the title to put in the title bar
3689 .  x, y - the screen coordinates of the upper left coordinate of the window
3690 .  m, n - the screen width and height in pixels
3691 -  howoften - if positive then determines the frequency of the plotting, if -1 then only at the final time
3692 
3693    Output Parameter:
3694 .  ctx - the context
3695 
3696    Options Database Key:
3697 +  -ts_monitor_lg_timestep - automatically sets line graph monitor
3698 .  -ts_monitor_lg_solution - monitor the solution (or certain values of the solution by calling TSMonitorLGSetDisplayVariables() or TSMonitorLGCtxSetDisplayVariables())
3699 .  -ts_monitor_lg_error -  monitor the error
3700 .  -ts_monitor_lg_ksp_iterations - monitor the number of KSP iterations needed for each timestep
3701 .  -ts_monitor_lg_snes_iterations - monitor the number of SNES iterations needed for each timestep
3702 -  -lg_use_markers <true,false> - mark the data points (at each time step) on the plot; default is true
3703 
3704    Notes:
3705    Use TSMonitorLGCtxDestroy() to destroy.
3706 
3707    One can provide a function that transforms the solution before plotting it with TSMonitorLGCtxSetTransform() or TSMonitorLGSetTransform()
3708 
3709    Many of the functions that control the monitoring have two forms: TSMonitorLGSet/GetXXXX() and TSMonitorLGCtxSet/GetXXXX() the first take a TS object as the
3710    first argument (if that TS object does not have a TSMonitorLGCtx associated with it the function call is ignored) and the second takes a TSMonitorLGCtx object
3711    as the first argument.
3712 
3713    One can control the names displayed for each solution or error variable with TSMonitorLGCtxSetVariableNames() or TSMonitorLGSetVariableNames()
3714 
3715 
3716    Level: intermediate
3717 
3718 .keywords: TS, monitor, line graph, residual
3719 
3720 .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError(), TSMonitorDefault(), VecView(),
3721            TSMonitorLGCtxCreate(), TSMonitorLGCtxSetVariableNames(), TSMonitorLGCtxGetVariableNames(),
3722            TSMonitorLGSetVariableNames(), TSMonitorLGGetVariableNames(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetDisplayVariables(),
3723            TSMonitorLGCtxSetTransform(), TSMonitorLGSetTransform(), TSMonitorLGError(), TSMonitorLGSNESIterations(), TSMonitorLGKSPIterations(),
3724            TSMonitorEnvelopeCtxCreate(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxDestroy(), TSMonitorEnvelop()
3725 
3726 @*/
3727 PetscErrorCode  TSMonitorLGCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtx *ctx)
3728 {
3729   PetscDraw      draw;
3730   PetscErrorCode ierr;
3731 
3732   PetscFunctionBegin;
3733   ierr = PetscNew(ctx);CHKERRQ(ierr);
3734   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&draw);CHKERRQ(ierr);
3735   ierr = PetscDrawSetFromOptions(draw);CHKERRQ(ierr);
3736   ierr = PetscDrawLGCreate(draw,1,&(*ctx)->lg);CHKERRQ(ierr);
3737   ierr = PetscDrawLGSetUseMarkers((*ctx)->lg,PETSC_TRUE);CHKERRQ(ierr);
3738   ierr = PetscDrawLGSetFromOptions((*ctx)->lg);CHKERRQ(ierr);
3739   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
3740   (*ctx)->howoften = howoften;
3741   PetscFunctionReturn(0);
3742 }
3743 
3744 #undef __FUNCT__
3745 #define __FUNCT__ "TSMonitorLGTimeStep"
3746 PetscErrorCode TSMonitorLGTimeStep(TS ts,PetscInt step,PetscReal ptime,Vec v,void *monctx)
3747 {
3748   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
3749   PetscReal      x   = ptime,y;
3750   PetscErrorCode ierr;
3751 
3752   PetscFunctionBegin;
3753   if (!step) {
3754     PetscDrawAxis axis;
3755     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
3756     ierr = PetscDrawAxisSetLabels(axis,"Timestep as function of time","Time","Time step");CHKERRQ(ierr);
3757     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
3758   }
3759   ierr = TSGetTimeStep(ts,&y);CHKERRQ(ierr);
3760   ierr = PetscDrawLGAddPoint(ctx->lg,&x,&y);CHKERRQ(ierr);
3761   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
3762     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
3763   }
3764   PetscFunctionReturn(0);
3765 }
3766 
3767 #undef __FUNCT__
3768 #define __FUNCT__ "TSMonitorLGCtxDestroy"
3769 /*@C
3770    TSMonitorLGCtxDestroy - Destroys a line graph context that was created
3771    with TSMonitorLGCtxCreate().
3772 
3773    Collective on TSMonitorLGCtx
3774 
3775    Input Parameter:
3776 .  ctx - the monitor context
3777 
3778    Level: intermediate
3779 
3780 .keywords: TS, monitor, line graph, destroy
3781 
3782 .seealso: TSMonitorLGCtxCreate(),  TSMonitorSet(), TSMonitorLGTimeStep();
3783 @*/
3784 PetscErrorCode  TSMonitorLGCtxDestroy(TSMonitorLGCtx *ctx)
3785 {
3786   PetscErrorCode ierr;
3787 
3788   PetscFunctionBegin;
3789   if ((*ctx)->transformdestroy) {
3790     ierr = ((*ctx)->transformdestroy)((*ctx)->transformctx);CHKERRQ(ierr);
3791   }
3792   ierr = PetscDrawLGDestroy(&(*ctx)->lg);CHKERRQ(ierr);
3793   ierr = PetscStrArrayDestroy(&(*ctx)->names);CHKERRQ(ierr);
3794   ierr = PetscStrArrayDestroy(&(*ctx)->displaynames);CHKERRQ(ierr);
3795   ierr = PetscFree((*ctx)->displayvariables);CHKERRQ(ierr);
3796   ierr = PetscFree((*ctx)->displayvalues);CHKERRQ(ierr);
3797   ierr = PetscFree(*ctx);CHKERRQ(ierr);
3798   PetscFunctionReturn(0);
3799 }
3800 
3801 #undef __FUNCT__
3802 #define __FUNCT__ "TSGetTime"
3803 /*@
3804    TSGetTime - Gets the time of the most recently completed step.
3805 
3806    Not Collective
3807 
3808    Input Parameter:
3809 .  ts - the TS context obtained from TSCreate()
3810 
3811    Output Parameter:
3812 .  t  - the current time
3813 
3814    Level: beginner
3815 
3816    Note:
3817    When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
3818    TSSetPreStage(), TSSetPostStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.
3819 
3820 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
3821 
3822 .keywords: TS, get, time
3823 @*/
3824 PetscErrorCode  TSGetTime(TS ts,PetscReal *t)
3825 {
3826   PetscFunctionBegin;
3827   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3828   PetscValidRealPointer(t,2);
3829   *t = ts->ptime;
3830   PetscFunctionReturn(0);
3831 }
3832 
3833 #undef __FUNCT__
3834 #define __FUNCT__ "TSGetPrevTime"
3835 /*@
3836    TSGetPrevTime - Gets the starting time of the previously completed step.
3837 
3838    Not Collective
3839 
3840    Input Parameter:
3841 .  ts - the TS context obtained from TSCreate()
3842 
3843    Output Parameter:
3844 .  t  - the previous time
3845 
3846    Level: beginner
3847 
3848 .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
3849 
3850 .keywords: TS, get, time
3851 @*/
3852 PetscErrorCode  TSGetPrevTime(TS ts,PetscReal *t)
3853 {
3854   PetscFunctionBegin;
3855   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3856   PetscValidRealPointer(t,2);
3857   *t = ts->ptime_prev;
3858   PetscFunctionReturn(0);
3859 }
3860 
3861 #undef __FUNCT__
3862 #define __FUNCT__ "TSSetTime"
3863 /*@
3864    TSSetTime - Allows one to reset the time.
3865 
3866    Logically Collective on TS
3867 
3868    Input Parameters:
3869 +  ts - the TS context obtained from TSCreate()
3870 -  time - the time
3871 
3872    Level: intermediate
3873 
3874 .seealso: TSGetTime(), TSSetDuration()
3875 
3876 .keywords: TS, set, time
3877 @*/
3878 PetscErrorCode  TSSetTime(TS ts, PetscReal t)
3879 {
3880   PetscFunctionBegin;
3881   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3882   PetscValidLogicalCollectiveReal(ts,t,2);
3883   ts->ptime = t;
3884   PetscFunctionReturn(0);
3885 }
3886 
3887 #undef __FUNCT__
3888 #define __FUNCT__ "TSSetOptionsPrefix"
3889 /*@C
3890    TSSetOptionsPrefix - Sets the prefix used for searching for all
3891    TS options in the database.
3892 
3893    Logically Collective on TS
3894 
3895    Input Parameter:
3896 +  ts     - The TS context
3897 -  prefix - The prefix to prepend to all option names
3898 
3899    Notes:
3900    A hyphen (-) must NOT be given at the beginning of the prefix name.
3901    The first character of all runtime options is AUTOMATICALLY the
3902    hyphen.
3903 
3904    Level: advanced
3905 
3906 .keywords: TS, set, options, prefix, database
3907 
3908 .seealso: TSSetFromOptions()
3909 
3910 @*/
3911 PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
3912 {
3913   PetscErrorCode ierr;
3914   SNES           snes;
3915 
3916   PetscFunctionBegin;
3917   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3918   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
3919   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3920   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
3921   PetscFunctionReturn(0);
3922 }
3923 
3924 
3925 #undef __FUNCT__
3926 #define __FUNCT__ "TSAppendOptionsPrefix"
3927 /*@C
3928    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
3929    TS options in the database.
3930 
3931    Logically Collective on TS
3932 
3933    Input Parameter:
3934 +  ts     - The TS context
3935 -  prefix - The prefix to prepend to all option names
3936 
3937    Notes:
3938    A hyphen (-) must NOT be given at the beginning of the prefix name.
3939    The first character of all runtime options is AUTOMATICALLY the
3940    hyphen.
3941 
3942    Level: advanced
3943 
3944 .keywords: TS, append, options, prefix, database
3945 
3946 .seealso: TSGetOptionsPrefix()
3947 
3948 @*/
3949 PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
3950 {
3951   PetscErrorCode ierr;
3952   SNES           snes;
3953 
3954   PetscFunctionBegin;
3955   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3956   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
3957   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3958   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
3959   PetscFunctionReturn(0);
3960 }
3961 
3962 #undef __FUNCT__
3963 #define __FUNCT__ "TSGetOptionsPrefix"
3964 /*@C
3965    TSGetOptionsPrefix - Sets the prefix used for searching for all
3966    TS options in the database.
3967 
3968    Not Collective
3969 
3970    Input Parameter:
3971 .  ts - The TS context
3972 
3973    Output Parameter:
3974 .  prefix - A pointer to the prefix string used
3975 
3976    Notes: On the fortran side, the user should pass in a string 'prifix' of
3977    sufficient length to hold the prefix.
3978 
3979    Level: intermediate
3980 
3981 .keywords: TS, get, options, prefix, database
3982 
3983 .seealso: TSAppendOptionsPrefix()
3984 @*/
3985 PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
3986 {
3987   PetscErrorCode ierr;
3988 
3989   PetscFunctionBegin;
3990   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3991   PetscValidPointer(prefix,2);
3992   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
3993   PetscFunctionReturn(0);
3994 }
3995 
3996 #undef __FUNCT__
3997 #define __FUNCT__ "TSGetRHSJacobian"
3998 /*@C
3999    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
4000 
4001    Not Collective, but parallel objects are returned if TS is parallel
4002 
4003    Input Parameter:
4004 .  ts  - The TS context obtained from TSCreate()
4005 
4006    Output Parameters:
4007 +  Amat - The (approximate) Jacobian J of G, where U_t = G(U,t)  (or NULL)
4008 .  Pmat - The matrix from which the preconditioner is constructed, usually the same as Amat  (or NULL)
4009 .  func - Function to compute the Jacobian of the RHS  (or NULL)
4010 -  ctx - User-defined context for Jacobian evaluation routine  (or NULL)
4011 
4012    Notes: You can pass in NULL for any return argument you do not need.
4013 
4014    Level: intermediate
4015 
4016 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
4017 
4018 .keywords: TS, timestep, get, matrix, Jacobian
4019 @*/
4020 PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *Amat,Mat *Pmat,TSRHSJacobian *func,void **ctx)
4021 {
4022   PetscErrorCode ierr;
4023   SNES           snes;
4024   DM             dm;
4025 
4026   PetscFunctionBegin;
4027   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
4028   ierr = SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);CHKERRQ(ierr);
4029   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
4030   ierr = DMTSGetRHSJacobian(dm,func,ctx);CHKERRQ(ierr);
4031   PetscFunctionReturn(0);
4032 }
4033 
4034 #undef __FUNCT__
4035 #define __FUNCT__ "TSGetIJacobian"
4036 /*@C
4037    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
4038 
4039    Not Collective, but parallel objects are returned if TS is parallel
4040 
4041    Input Parameter:
4042 .  ts  - The TS context obtained from TSCreate()
4043 
4044    Output Parameters:
4045 +  Amat  - The (approximate) Jacobian of F(t,U,U_t)
4046 .  Pmat - The matrix from which the preconditioner is constructed, often the same as Amat
4047 .  f   - The function to compute the matrices
4048 - ctx - User-defined context for Jacobian evaluation routine
4049 
4050    Notes: You can pass in NULL for any return argument you do not need.
4051 
4052    Level: advanced
4053 
4054 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
4055 
4056 .keywords: TS, timestep, get, matrix, Jacobian
4057 @*/
4058 PetscErrorCode  TSGetIJacobian(TS ts,Mat *Amat,Mat *Pmat,TSIJacobian *f,void **ctx)
4059 {
4060   PetscErrorCode ierr;
4061   SNES           snes;
4062   DM             dm;
4063 
4064   PetscFunctionBegin;
4065   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
4066   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
4067   ierr = SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);CHKERRQ(ierr);
4068   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
4069   ierr = DMTSGetIJacobian(dm,f,ctx);CHKERRQ(ierr);
4070   PetscFunctionReturn(0);
4071 }
4072 
4073 
4074 #undef __FUNCT__
4075 #define __FUNCT__ "TSMonitorDrawSolution"
4076 /*@C
4077    TSMonitorDrawSolution - Monitors progress of the TS solvers by calling
4078    VecView() for the solution at each timestep
4079 
4080    Collective on TS
4081 
4082    Input Parameters:
4083 +  ts - the TS context
4084 .  step - current time-step
4085 .  ptime - current time
4086 -  dummy - either a viewer or NULL
4087 
4088    Options Database:
4089 .   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
4090 
4091    Notes: the initial solution and current solution are not display with a common axis scaling so generally the option -ts_monitor_draw_solution_initial
4092        will look bad
4093 
4094    Level: intermediate
4095 
4096 .keywords: TS,  vector, monitor, view
4097 
4098 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4099 @*/
4100 PetscErrorCode  TSMonitorDrawSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4101 {
4102   PetscErrorCode   ierr;
4103   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
4104   PetscDraw        draw;
4105 
4106   PetscFunctionBegin;
4107   if (!step && ictx->showinitial) {
4108     if (!ictx->initialsolution) {
4109       ierr = VecDuplicate(u,&ictx->initialsolution);CHKERRQ(ierr);
4110     }
4111     ierr = VecCopy(u,ictx->initialsolution);CHKERRQ(ierr);
4112   }
4113   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
4114 
4115   if (ictx->showinitial) {
4116     PetscReal pause;
4117     ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr);
4118     ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr);
4119     ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr);
4120     ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr);
4121     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr);
4122   }
4123   ierr = VecView(u,ictx->viewer);CHKERRQ(ierr);
4124   if (ictx->showtimestepandtime) {
4125     PetscReal xl,yl,xr,yr,h;
4126     char      time[32];
4127 
4128     ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
4129     ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
4130     ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
4131     h    = yl + .95*(yr - yl);
4132     ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
4133     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
4134   }
4135 
4136   if (ictx->showinitial) {
4137     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr);
4138   }
4139   PetscFunctionReturn(0);
4140 }
4141 
4142 #undef __FUNCT__
4143 #define __FUNCT__ "TSAdjointMonitorDrawSensi"
4144 /*@C
4145    TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling
4146    VecView() for the sensitivities to initial states at each timestep
4147 
4148    Collective on TS
4149 
4150    Input Parameters:
4151 +  ts - the TS context
4152 .  step - current time-step
4153 .  ptime - current time
4154 .  u - current state
4155 .  numcost - number of cost functions
4156 .  lambda - sensitivities to initial conditions
4157 .  mu - sensitivities to parameters
4158 -  dummy - either a viewer or NULL
4159 
4160    Level: intermediate
4161 
4162 .keywords: TS,  vector, adjoint, monitor, view
4163 
4164 .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView()
4165 @*/
4166 PetscErrorCode  TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
4167 {
4168   PetscErrorCode   ierr;
4169   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
4170   PetscDraw        draw;
4171   PetscReal        xl,yl,xr,yr,h;
4172   char             time[32];
4173 
4174   PetscFunctionBegin;
4175   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
4176 
4177   ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr);
4178   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
4179   ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
4180   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
4181   h    = yl + .95*(yr - yl);
4182   ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
4183   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
4184 
4185   PetscFunctionReturn(0);
4186 }
4187 
4188 #undef __FUNCT__
4189 #define __FUNCT__ "TSMonitorDrawSolutionPhase"
4190 /*@C
4191    TSMonitorDrawSolutionPhase - Monitors progress of the TS solvers by plotting the solution as a phase diagram
4192 
4193    Collective on TS
4194 
4195    Input Parameters:
4196 +  ts - the TS context
4197 .  step - current time-step
4198 .  ptime - current time
4199 -  dummy - either a viewer or NULL
4200 
4201    Level: intermediate
4202 
4203 .keywords: TS,  vector, monitor, view
4204 
4205 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4206 @*/
4207 PetscErrorCode  TSMonitorDrawSolutionPhase(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4208 {
4209   PetscErrorCode    ierr;
4210   TSMonitorDrawCtx  ictx = (TSMonitorDrawCtx)dummy;
4211   PetscDraw         draw;
4212   MPI_Comm          comm;
4213   PetscInt          n;
4214   PetscMPIInt       size;
4215   PetscReal         xl,yl,xr,yr,h;
4216   char              time[32];
4217   const PetscScalar *U;
4218 
4219   PetscFunctionBegin;
4220   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
4221   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
4222   if (size != 1) SETERRQ(comm,PETSC_ERR_SUP,"Only allowed for sequential runs");
4223   ierr = VecGetSize(u,&n);CHKERRQ(ierr);
4224   if (n != 2) SETERRQ(comm,PETSC_ERR_SUP,"Only for ODEs with two unknowns");
4225 
4226   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
4227 
4228   ierr = VecGetArrayRead(u,&U);CHKERRQ(ierr);
4229   ierr = PetscDrawAxisGetLimits(ictx->axis,&xl,&xr,&yl,&yr);CHKERRQ(ierr);
4230   if ((PetscRealPart(U[0]) < xl) || (PetscRealPart(U[1]) < yl) || (PetscRealPart(U[0]) > xr) || (PetscRealPart(U[1]) > yr)) {
4231       ierr = VecRestoreArrayRead(u,&U);CHKERRQ(ierr);
4232       PetscFunctionReturn(0);
4233   }
4234   if (!step) ictx->color++;
4235   ierr = PetscDrawPoint(draw,PetscRealPart(U[0]),PetscRealPart(U[1]),ictx->color);CHKERRQ(ierr);
4236   ierr = VecRestoreArrayRead(u,&U);CHKERRQ(ierr);
4237 
4238   if (ictx->showtimestepandtime) {
4239     ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
4240     ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
4241     h    = yl + .95*(yr - yl);
4242     ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
4243   }
4244   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
4245   PetscFunctionReturn(0);
4246 }
4247 
4248 
4249 #undef __FUNCT__
4250 #define __FUNCT__ "TSMonitorDrawCtxDestroy"
4251 /*@C
4252    TSMonitorDrawCtxDestroy - Destroys the monitor context for TSMonitorDrawSolution()
4253 
4254    Collective on TS
4255 
4256    Input Parameters:
4257 .    ctx - the monitor context
4258 
4259    Level: intermediate
4260 
4261 .keywords: TS,  vector, monitor, view
4262 
4263 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawSolution(), TSMonitorDrawError()
4264 @*/
4265 PetscErrorCode  TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *ictx)
4266 {
4267   PetscErrorCode ierr;
4268 
4269   PetscFunctionBegin;
4270   ierr = PetscDrawAxisDestroy(&(*ictx)->axis);CHKERRQ(ierr);
4271   ierr = PetscViewerDestroy(&(*ictx)->viewer);CHKERRQ(ierr);
4272   ierr = VecDestroy(&(*ictx)->initialsolution);CHKERRQ(ierr);
4273   ierr = PetscFree(*ictx);CHKERRQ(ierr);
4274   PetscFunctionReturn(0);
4275 }
4276 
4277 #undef __FUNCT__
4278 #define __FUNCT__ "TSMonitorDrawCtxCreate"
4279 /*@C
4280    TSMonitorDrawCtxCreate - Creates the monitor context for TSMonitorDrawCtx
4281 
4282    Collective on TS
4283 
4284    Input Parameter:
4285 .    ts - time-step context
4286 
4287    Output Patameter:
4288 .    ctx - the monitor context
4289 
4290    Options Database:
4291 .   -ts_monitor_draw_solution_initial - show initial solution as well as current solution
4292 
4293    Level: intermediate
4294 
4295 .keywords: TS,  vector, monitor, view
4296 
4297 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorDrawCtx()
4298 @*/
4299 PetscErrorCode  TSMonitorDrawCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorDrawCtx *ctx)
4300 {
4301   PetscErrorCode   ierr;
4302 
4303   PetscFunctionBegin;
4304   ierr = PetscNew(ctx);CHKERRQ(ierr);
4305   ierr = PetscViewerDrawOpen(comm,host,label,x,y,m,n,&(*ctx)->viewer);CHKERRQ(ierr);
4306   ierr = PetscViewerSetFromOptions((*ctx)->viewer);CHKERRQ(ierr);
4307 
4308   (*ctx)->howoften    = howoften;
4309   (*ctx)->showinitial = PETSC_FALSE;
4310   ierr = PetscOptionsGetBool(NULL,NULL,"-ts_monitor_draw_solution_initial",&(*ctx)->showinitial,NULL);CHKERRQ(ierr);
4311 
4312   (*ctx)->showtimestepandtime = PETSC_FALSE;
4313   ierr = PetscOptionsGetBool(NULL,NULL,"-ts_monitor_draw_solution_show_time",&(*ctx)->showtimestepandtime,NULL);CHKERRQ(ierr);
4314   (*ctx)->color = PETSC_DRAW_WHITE;
4315   PetscFunctionReturn(0);
4316 }
4317 
4318 #undef __FUNCT__
4319 #define __FUNCT__ "TSMonitorDrawError"
4320 /*@C
4321    TSMonitorDrawError - Monitors progress of the TS solvers by calling
4322    VecView() for the error at each timestep
4323 
4324    Collective on TS
4325 
4326    Input Parameters:
4327 +  ts - the TS context
4328 .  step - current time-step
4329 .  ptime - current time
4330 -  dummy - either a viewer or NULL
4331 
4332    Level: intermediate
4333 
4334 .keywords: TS,  vector, monitor, view
4335 
4336 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4337 @*/
4338 PetscErrorCode  TSMonitorDrawError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
4339 {
4340   PetscErrorCode   ierr;
4341   TSMonitorDrawCtx ctx    = (TSMonitorDrawCtx)dummy;
4342   PetscViewer      viewer = ctx->viewer;
4343   Vec              work;
4344 
4345   PetscFunctionBegin;
4346   if (!(((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
4347   ierr = VecDuplicate(u,&work);CHKERRQ(ierr);
4348   ierr = TSComputeSolutionFunction(ts,ptime,work);CHKERRQ(ierr);
4349   ierr = VecAXPY(work,-1.0,u);CHKERRQ(ierr);
4350   ierr = VecView(work,viewer);CHKERRQ(ierr);
4351   ierr = VecDestroy(&work);CHKERRQ(ierr);
4352   PetscFunctionReturn(0);
4353 }
4354 
4355 #include <petsc/private/dmimpl.h>
4356 #undef __FUNCT__
4357 #define __FUNCT__ "TSSetDM"
4358 /*@
4359    TSSetDM - Sets the DM that may be used by some preconditioners
4360 
4361    Logically Collective on TS and DM
4362 
4363    Input Parameters:
4364 +  ts - the preconditioner context
4365 -  dm - the dm
4366 
4367    Level: intermediate
4368 
4369 
4370 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
4371 @*/
4372 PetscErrorCode  TSSetDM(TS ts,DM dm)
4373 {
4374   PetscErrorCode ierr;
4375   SNES           snes;
4376   DMTS           tsdm;
4377 
4378   PetscFunctionBegin;
4379   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4380   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
4381   if (ts->dm) {               /* Move the DMTS context over to the new DM unless the new DM already has one */
4382     if (ts->dm->dmts && !dm->dmts) {
4383       ierr = DMCopyDMTS(ts->dm,dm);CHKERRQ(ierr);
4384       ierr = DMGetDMTS(ts->dm,&tsdm);CHKERRQ(ierr);
4385       if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
4386         tsdm->originaldm = dm;
4387       }
4388     }
4389     ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
4390   }
4391   ts->dm = dm;
4392 
4393   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
4394   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
4395   PetscFunctionReturn(0);
4396 }
4397 
4398 #undef __FUNCT__
4399 #define __FUNCT__ "TSGetDM"
4400 /*@
4401    TSGetDM - Gets the DM that may be used by some preconditioners
4402 
4403    Not Collective
4404 
4405    Input Parameter:
4406 . ts - the preconditioner context
4407 
4408    Output Parameter:
4409 .  dm - the dm
4410 
4411    Level: intermediate
4412 
4413 
4414 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
4415 @*/
4416 PetscErrorCode  TSGetDM(TS ts,DM *dm)
4417 {
4418   PetscErrorCode ierr;
4419 
4420   PetscFunctionBegin;
4421   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4422   if (!ts->dm) {
4423     ierr = DMShellCreate(PetscObjectComm((PetscObject)ts),&ts->dm);CHKERRQ(ierr);
4424     if (ts->snes) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
4425   }
4426   *dm = ts->dm;
4427   PetscFunctionReturn(0);
4428 }
4429 
4430 #undef __FUNCT__
4431 #define __FUNCT__ "SNESTSFormFunction"
4432 /*@
4433    SNESTSFormFunction - Function to evaluate nonlinear residual
4434 
4435    Logically Collective on SNES
4436 
4437    Input Parameter:
4438 + snes - nonlinear solver
4439 . U - the current state at which to evaluate the residual
4440 - ctx - user context, must be a TS
4441 
4442    Output Parameter:
4443 . F - the nonlinear residual
4444 
4445    Notes:
4446    This function is not normally called by users and is automatically registered with the SNES used by TS.
4447    It is most frequently passed to MatFDColoringSetFunction().
4448 
4449    Level: advanced
4450 
4451 .seealso: SNESSetFunction(), MatFDColoringSetFunction()
4452 @*/
4453 PetscErrorCode  SNESTSFormFunction(SNES snes,Vec U,Vec F,void *ctx)
4454 {
4455   TS             ts = (TS)ctx;
4456   PetscErrorCode ierr;
4457 
4458   PetscFunctionBegin;
4459   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4460   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
4461   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
4462   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
4463   ierr = (ts->ops->snesfunction)(snes,U,F,ts);CHKERRQ(ierr);
4464   PetscFunctionReturn(0);
4465 }
4466 
4467 #undef __FUNCT__
4468 #define __FUNCT__ "SNESTSFormJacobian"
4469 /*@
4470    SNESTSFormJacobian - Function to evaluate the Jacobian
4471 
4472    Collective on SNES
4473 
4474    Input Parameter:
4475 + snes - nonlinear solver
4476 . U - the current state at which to evaluate the residual
4477 - ctx - user context, must be a TS
4478 
4479    Output Parameter:
4480 + A - the Jacobian
4481 . B - the preconditioning matrix (may be the same as A)
4482 - flag - indicates any structure change in the matrix
4483 
4484    Notes:
4485    This function is not normally called by users and is automatically registered with the SNES used by TS.
4486 
4487    Level: developer
4488 
4489 .seealso: SNESSetJacobian()
4490 @*/
4491 PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec U,Mat A,Mat B,void *ctx)
4492 {
4493   TS             ts = (TS)ctx;
4494   PetscErrorCode ierr;
4495 
4496   PetscFunctionBegin;
4497   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4498   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
4499   PetscValidPointer(A,3);
4500   PetscValidHeaderSpecific(A,MAT_CLASSID,3);
4501   PetscValidPointer(B,4);
4502   PetscValidHeaderSpecific(B,MAT_CLASSID,4);
4503   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
4504   ierr = (ts->ops->snesjacobian)(snes,U,A,B,ts);CHKERRQ(ierr);
4505   PetscFunctionReturn(0);
4506 }
4507 
4508 #undef __FUNCT__
4509 #define __FUNCT__ "TSComputeRHSFunctionLinear"
4510 /*@C
4511    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems Udot = A U only
4512 
4513    Collective on TS
4514 
4515    Input Arguments:
4516 +  ts - time stepping context
4517 .  t - time at which to evaluate
4518 .  U - state at which to evaluate
4519 -  ctx - context
4520 
4521    Output Arguments:
4522 .  F - right hand side
4523 
4524    Level: intermediate
4525 
4526    Notes:
4527    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
4528    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
4529 
4530 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
4531 @*/
4532 PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
4533 {
4534   PetscErrorCode ierr;
4535   Mat            Arhs,Brhs;
4536 
4537   PetscFunctionBegin;
4538   ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
4539   ierr = TSComputeRHSJacobian(ts,t,U,Arhs,Brhs);CHKERRQ(ierr);
4540   ierr = MatMult(Arhs,U,F);CHKERRQ(ierr);
4541   PetscFunctionReturn(0);
4542 }
4543 
4544 #undef __FUNCT__
4545 #define __FUNCT__ "TSComputeRHSJacobianConstant"
4546 /*@C
4547    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
4548 
4549    Collective on TS
4550 
4551    Input Arguments:
4552 +  ts - time stepping context
4553 .  t - time at which to evaluate
4554 .  U - state at which to evaluate
4555 -  ctx - context
4556 
4557    Output Arguments:
4558 +  A - pointer to operator
4559 .  B - pointer to preconditioning matrix
4560 -  flg - matrix structure flag
4561 
4562    Level: intermediate
4563 
4564    Notes:
4565    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
4566 
4567 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
4568 @*/
4569 PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
4570 {
4571   PetscFunctionBegin;
4572   PetscFunctionReturn(0);
4573 }
4574 
4575 #undef __FUNCT__
4576 #define __FUNCT__ "TSComputeIFunctionLinear"
4577 /*@C
4578    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
4579 
4580    Collective on TS
4581 
4582    Input Arguments:
4583 +  ts - time stepping context
4584 .  t - time at which to evaluate
4585 .  U - state at which to evaluate
4586 .  Udot - time derivative of state vector
4587 -  ctx - context
4588 
4589    Output Arguments:
4590 .  F - left hand side
4591 
4592    Level: intermediate
4593 
4594    Notes:
4595    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
4596    user is required to write their own TSComputeIFunction.
4597    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
4598    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
4599 
4600    Note that using this function is NOT equivalent to using TSComputeRHSFunctionLinear() since that solves Udot = A U
4601 
4602 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant(), TSComputeRHSFunctionLinear()
4603 @*/
4604 PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
4605 {
4606   PetscErrorCode ierr;
4607   Mat            A,B;
4608 
4609   PetscFunctionBegin;
4610   ierr = TSGetIJacobian(ts,&A,&B,NULL,NULL);CHKERRQ(ierr);
4611   ierr = TSComputeIJacobian(ts,t,U,Udot,1.0,A,B,PETSC_TRUE);CHKERRQ(ierr);
4612   ierr = MatMult(A,Udot,F);CHKERRQ(ierr);
4613   PetscFunctionReturn(0);
4614 }
4615 
4616 #undef __FUNCT__
4617 #define __FUNCT__ "TSComputeIJacobianConstant"
4618 /*@C
4619    TSComputeIJacobianConstant - Reuses a time-independent for a semi-implicit DAE or ODE
4620 
4621    Collective on TS
4622 
4623    Input Arguments:
4624 +  ts - time stepping context
4625 .  t - time at which to evaluate
4626 .  U - state at which to evaluate
4627 .  Udot - time derivative of state vector
4628 .  shift - shift to apply
4629 -  ctx - context
4630 
4631    Output Arguments:
4632 +  A - pointer to operator
4633 .  B - pointer to preconditioning matrix
4634 -  flg - matrix structure flag
4635 
4636    Level: advanced
4637 
4638    Notes:
4639    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
4640 
4641    It is only appropriate for problems of the form
4642 
4643 $     M Udot = F(U,t)
4644 
4645   where M is constant and F is non-stiff.  The user must pass M to TSSetIJacobian().  The current implementation only
4646   works with IMEX time integration methods such as TSROSW and TSARKIMEX, since there is no support for de-constructing
4647   an implicit operator of the form
4648 
4649 $    shift*M + J
4650 
4651   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
4652   a copy of M or reassemble it when requested.
4653 
4654 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
4655 @*/
4656 PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat A,Mat B,void *ctx)
4657 {
4658   PetscErrorCode ierr;
4659 
4660   PetscFunctionBegin;
4661   ierr = MatScale(A, shift / ts->ijacobian.shift);CHKERRQ(ierr);
4662   ts->ijacobian.shift = shift;
4663   PetscFunctionReturn(0);
4664 }
4665 
4666 #undef __FUNCT__
4667 #define __FUNCT__ "TSGetEquationType"
4668 /*@
4669    TSGetEquationType - Gets the type of the equation that TS is solving.
4670 
4671    Not Collective
4672 
4673    Input Parameter:
4674 .  ts - the TS context
4675 
4676    Output Parameter:
4677 .  equation_type - see TSEquationType
4678 
4679    Level: beginner
4680 
4681 .keywords: TS, equation type
4682 
4683 .seealso: TSSetEquationType(), TSEquationType
4684 @*/
4685 PetscErrorCode  TSGetEquationType(TS ts,TSEquationType *equation_type)
4686 {
4687   PetscFunctionBegin;
4688   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4689   PetscValidPointer(equation_type,2);
4690   *equation_type = ts->equation_type;
4691   PetscFunctionReturn(0);
4692 }
4693 
4694 #undef __FUNCT__
4695 #define __FUNCT__ "TSSetEquationType"
4696 /*@
4697    TSSetEquationType - Sets the type of the equation that TS is solving.
4698 
4699    Not Collective
4700 
4701    Input Parameter:
4702 +  ts - the TS context
4703 -  equation_type - see TSEquationType
4704 
4705    Level: advanced
4706 
4707 .keywords: TS, equation type
4708 
4709 .seealso: TSGetEquationType(), TSEquationType
4710 @*/
4711 PetscErrorCode  TSSetEquationType(TS ts,TSEquationType equation_type)
4712 {
4713   PetscFunctionBegin;
4714   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4715   ts->equation_type = equation_type;
4716   PetscFunctionReturn(0);
4717 }
4718 
4719 #undef __FUNCT__
4720 #define __FUNCT__ "TSGetConvergedReason"
4721 /*@
4722    TSGetConvergedReason - Gets the reason the TS iteration was stopped.
4723 
4724    Not Collective
4725 
4726    Input Parameter:
4727 .  ts - the TS context
4728 
4729    Output Parameter:
4730 .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
4731             manual pages for the individual convergence tests for complete lists
4732 
4733    Level: beginner
4734 
4735    Notes:
4736    Can only be called after the call to TSSolve() is complete.
4737 
4738 .keywords: TS, nonlinear, set, convergence, test
4739 
4740 .seealso: TSSetConvergenceTest(), TSConvergedReason
4741 @*/
4742 PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
4743 {
4744   PetscFunctionBegin;
4745   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4746   PetscValidPointer(reason,2);
4747   *reason = ts->reason;
4748   PetscFunctionReturn(0);
4749 }
4750 
4751 #undef __FUNCT__
4752 #define __FUNCT__ "TSSetConvergedReason"
4753 /*@
4754    TSSetConvergedReason - Sets the reason for handling the convergence of TSSolve.
4755 
4756    Not Collective
4757 
4758    Input Parameter:
4759 +  ts - the TS context
4760 .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
4761             manual pages for the individual convergence tests for complete lists
4762 
4763    Level: advanced
4764 
4765    Notes:
4766    Can only be called during TSSolve() is active.
4767 
4768 .keywords: TS, nonlinear, set, convergence, test
4769 
4770 .seealso: TSConvergedReason
4771 @*/
4772 PetscErrorCode  TSSetConvergedReason(TS ts,TSConvergedReason reason)
4773 {
4774   PetscFunctionBegin;
4775   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4776   ts->reason = reason;
4777   PetscFunctionReturn(0);
4778 }
4779 
4780 #undef __FUNCT__
4781 #define __FUNCT__ "TSGetSolveTime"
4782 /*@
4783    TSGetSolveTime - Gets the time after a call to TSSolve()
4784 
4785    Not Collective
4786 
4787    Input Parameter:
4788 .  ts - the TS context
4789 
4790    Output Parameter:
4791 .  ftime - the final time. This time should correspond to the final time set with TSSetDuration()
4792 
4793    Level: beginner
4794 
4795    Notes:
4796    Can only be called after the call to TSSolve() is complete.
4797 
4798 .keywords: TS, nonlinear, set, convergence, test
4799 
4800 .seealso: TSSetConvergenceTest(), TSConvergedReason
4801 @*/
4802 PetscErrorCode  TSGetSolveTime(TS ts,PetscReal *ftime)
4803 {
4804   PetscFunctionBegin;
4805   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4806   PetscValidPointer(ftime,2);
4807   *ftime = ts->solvetime;
4808   PetscFunctionReturn(0);
4809 }
4810 
4811 #undef __FUNCT__
4812 #define __FUNCT__ "TSGetTotalSteps"
4813 /*@
4814    TSGetTotalSteps - Gets the total number of steps done since the last call to TSSetUp() or TSCreate()
4815 
4816    Not Collective
4817 
4818    Input Parameter:
4819 .  ts - the TS context
4820 
4821    Output Parameter:
4822 .  steps - the number of steps
4823 
4824    Level: beginner
4825 
4826    Notes:
4827    Includes the number of steps for all calls to TSSolve() since TSSetUp() was called
4828 
4829 .keywords: TS, nonlinear, set, convergence, test
4830 
4831 .seealso: TSSetConvergenceTest(), TSConvergedReason
4832 @*/
4833 PetscErrorCode  TSGetTotalSteps(TS ts,PetscInt *steps)
4834 {
4835   PetscFunctionBegin;
4836   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4837   PetscValidPointer(steps,2);
4838   *steps = ts->total_steps;
4839   PetscFunctionReturn(0);
4840 }
4841 
4842 #undef __FUNCT__
4843 #define __FUNCT__ "TSGetSNESIterations"
4844 /*@
4845    TSGetSNESIterations - Gets the total number of nonlinear iterations
4846    used by the time integrator.
4847 
4848    Not Collective
4849 
4850    Input Parameter:
4851 .  ts - TS context
4852 
4853    Output Parameter:
4854 .  nits - number of nonlinear iterations
4855 
4856    Notes:
4857    This counter is reset to zero for each successive call to TSSolve().
4858 
4859    Level: intermediate
4860 
4861 .keywords: TS, get, number, nonlinear, iterations
4862 
4863 .seealso:  TSGetKSPIterations()
4864 @*/
4865 PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
4866 {
4867   PetscFunctionBegin;
4868   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4869   PetscValidIntPointer(nits,2);
4870   *nits = ts->snes_its;
4871   PetscFunctionReturn(0);
4872 }
4873 
4874 #undef __FUNCT__
4875 #define __FUNCT__ "TSGetKSPIterations"
4876 /*@
4877    TSGetKSPIterations - Gets the total number of linear iterations
4878    used by the time integrator.
4879 
4880    Not Collective
4881 
4882    Input Parameter:
4883 .  ts - TS context
4884 
4885    Output Parameter:
4886 .  lits - number of linear iterations
4887 
4888    Notes:
4889    This counter is reset to zero for each successive call to TSSolve().
4890 
4891    Level: intermediate
4892 
4893 .keywords: TS, get, number, linear, iterations
4894 
4895 .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
4896 @*/
4897 PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
4898 {
4899   PetscFunctionBegin;
4900   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4901   PetscValidIntPointer(lits,2);
4902   *lits = ts->ksp_its;
4903   PetscFunctionReturn(0);
4904 }
4905 
4906 #undef __FUNCT__
4907 #define __FUNCT__ "TSGetStepRejections"
4908 /*@
4909    TSGetStepRejections - Gets the total number of rejected steps.
4910 
4911    Not Collective
4912 
4913    Input Parameter:
4914 .  ts - TS context
4915 
4916    Output Parameter:
4917 .  rejects - number of steps rejected
4918 
4919    Notes:
4920    This counter is reset to zero for each successive call to TSSolve().
4921 
4922    Level: intermediate
4923 
4924 .keywords: TS, get, number
4925 
4926 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
4927 @*/
4928 PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
4929 {
4930   PetscFunctionBegin;
4931   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4932   PetscValidIntPointer(rejects,2);
4933   *rejects = ts->reject;
4934   PetscFunctionReturn(0);
4935 }
4936 
4937 #undef __FUNCT__
4938 #define __FUNCT__ "TSGetSNESFailures"
4939 /*@
4940    TSGetSNESFailures - Gets the total number of failed SNES solves
4941 
4942    Not Collective
4943 
4944    Input Parameter:
4945 .  ts - TS context
4946 
4947    Output Parameter:
4948 .  fails - number of failed nonlinear solves
4949 
4950    Notes:
4951    This counter is reset to zero for each successive call to TSSolve().
4952 
4953    Level: intermediate
4954 
4955 .keywords: TS, get, number
4956 
4957 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
4958 @*/
4959 PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
4960 {
4961   PetscFunctionBegin;
4962   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4963   PetscValidIntPointer(fails,2);
4964   *fails = ts->num_snes_failures;
4965   PetscFunctionReturn(0);
4966 }
4967 
4968 #undef __FUNCT__
4969 #define __FUNCT__ "TSSetMaxStepRejections"
4970 /*@
4971    TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails
4972 
4973    Not Collective
4974 
4975    Input Parameter:
4976 +  ts - TS context
4977 -  rejects - maximum number of rejected steps, pass -1 for unlimited
4978 
4979    Notes:
4980    The counter is reset to zero for each step
4981 
4982    Options Database Key:
4983  .  -ts_max_reject - Maximum number of step rejections before a step fails
4984 
4985    Level: intermediate
4986 
4987 .keywords: TS, set, maximum, number
4988 
4989 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
4990 @*/
4991 PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
4992 {
4993   PetscFunctionBegin;
4994   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4995   ts->max_reject = rejects;
4996   PetscFunctionReturn(0);
4997 }
4998 
4999 #undef __FUNCT__
5000 #define __FUNCT__ "TSSetMaxSNESFailures"
5001 /*@
5002    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves
5003 
5004    Not Collective
5005 
5006    Input Parameter:
5007 +  ts - TS context
5008 -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited
5009 
5010    Notes:
5011    The counter is reset to zero for each successive call to TSSolve().
5012 
5013    Options Database Key:
5014  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures
5015 
5016    Level: intermediate
5017 
5018 .keywords: TS, set, maximum, number
5019 
5020 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
5021 @*/
5022 PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
5023 {
5024   PetscFunctionBegin;
5025   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5026   ts->max_snes_failures = fails;
5027   PetscFunctionReturn(0);
5028 }
5029 
5030 #undef __FUNCT__
5031 #define __FUNCT__ "TSSetErrorIfStepFails"
5032 /*@
5033    TSSetErrorIfStepFails - Error if no step succeeds
5034 
5035    Not Collective
5036 
5037    Input Parameter:
5038 +  ts - TS context
5039 -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure
5040 
5041    Options Database Key:
5042  .  -ts_error_if_step_fails - Error if no step succeeds
5043 
5044    Level: intermediate
5045 
5046 .keywords: TS, set, error
5047 
5048 .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
5049 @*/
5050 PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
5051 {
5052   PetscFunctionBegin;
5053   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5054   ts->errorifstepfailed = err;
5055   PetscFunctionReturn(0);
5056 }
5057 
5058 #undef __FUNCT__
5059 #define __FUNCT__ "TSMonitorSolution"
5060 /*@C
5061    TSMonitorSolution - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file or a PetscDraw object
5062 
5063    Collective on TS
5064 
5065    Input Parameters:
5066 +  ts - the TS context
5067 .  step - current time-step
5068 .  ptime - current time
5069 .  u - current state
5070 -  viewer - binary viewer
5071 
5072    Level: intermediate
5073 
5074 .keywords: TS,  vector, monitor, view
5075 
5076 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
5077 @*/
5078 PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *viewer)
5079 {
5080   PetscErrorCode ierr;
5081   PetscViewer    v = (PetscViewer)viewer;
5082 
5083   PetscFunctionBegin;
5084   ierr = VecView(u,v);CHKERRQ(ierr);
5085   PetscFunctionReturn(0);
5086 }
5087 
5088 #undef __FUNCT__
5089 #define __FUNCT__ "TSMonitorSolutionVTK"
5090 /*@C
5091    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.
5092 
5093    Collective on TS
5094 
5095    Input Parameters:
5096 +  ts - the TS context
5097 .  step - current time-step
5098 .  ptime - current time
5099 .  u - current state
5100 -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
5101 
5102    Level: intermediate
5103 
5104    Notes:
5105    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.
5106    These are named according to the file name template.
5107 
5108    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().
5109 
5110 .keywords: TS,  vector, monitor, view
5111 
5112 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
5113 @*/
5114 PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec u,void *filenametemplate)
5115 {
5116   PetscErrorCode ierr;
5117   char           filename[PETSC_MAX_PATH_LEN];
5118   PetscViewer    viewer;
5119 
5120   PetscFunctionBegin;
5121   ierr = PetscSNPrintf(filename,sizeof(filename),(const char*)filenametemplate,step);CHKERRQ(ierr);
5122   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)ts),filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
5123   ierr = VecView(u,viewer);CHKERRQ(ierr);
5124   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
5125   PetscFunctionReturn(0);
5126 }
5127 
5128 #undef __FUNCT__
5129 #define __FUNCT__ "TSMonitorSolutionVTKDestroy"
5130 /*@C
5131    TSMonitorSolutionVTKDestroy - Destroy context for monitoring
5132 
5133    Collective on TS
5134 
5135    Input Parameters:
5136 .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
5137 
5138    Level: intermediate
5139 
5140    Note:
5141    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().
5142 
5143 .keywords: TS,  vector, monitor, view
5144 
5145 .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
5146 @*/
5147 PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
5148 {
5149   PetscErrorCode ierr;
5150 
5151   PetscFunctionBegin;
5152   ierr = PetscFree(*(char**)filenametemplate);CHKERRQ(ierr);
5153   PetscFunctionReturn(0);
5154 }
5155 
5156 #undef __FUNCT__
5157 #define __FUNCT__ "TSGetAdapt"
5158 /*@
5159    TSGetAdapt - Get the adaptive controller context for the current method
5160 
5161    Collective on TS if controller has not been created yet
5162 
5163    Input Arguments:
5164 .  ts - time stepping context
5165 
5166    Output Arguments:
5167 .  adapt - adaptive controller
5168 
5169    Level: intermediate
5170 
5171 .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
5172 @*/
5173 PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
5174 {
5175   PetscErrorCode ierr;
5176 
5177   PetscFunctionBegin;
5178   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5179   PetscValidPointer(adapt,2);
5180   if (!ts->adapt) {
5181     ierr = TSAdaptCreate(PetscObjectComm((PetscObject)ts),&ts->adapt);CHKERRQ(ierr);
5182     ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->adapt);CHKERRQ(ierr);
5183     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);CHKERRQ(ierr);
5184   }
5185   *adapt = ts->adapt;
5186   PetscFunctionReturn(0);
5187 }
5188 
5189 #undef __FUNCT__
5190 #define __FUNCT__ "TSSetTolerances"
5191 /*@
5192    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller
5193 
5194    Logically Collective
5195 
5196    Input Arguments:
5197 +  ts - time integration context
5198 .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
5199 .  vatol - vector of absolute tolerances or NULL, used in preference to atol if present
5200 .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
5201 -  vrtol - vector of relative tolerances or NULL, used in preference to atol if present
5202 
5203    Options Database keys:
5204 +  -ts_rtol <rtol> - relative tolerance for local truncation error
5205 -  -ts_atol <atol> Absolute tolerance for local truncation error
5206 
5207    Notes:
5208    With PETSc's implicit schemes for DAE problems, the calculation of the local truncation error
5209    (LTE) includes both the differential and the algebraic variables. If one wants the LTE to be
5210    computed only for the differential or the algebraic part then this can be done using the vector of
5211    tolerances vatol. For example, by setting the tolerance vector with the desired tolerance for the
5212    differential part and infinity for the algebraic part, the LTE calculation will include only the
5213    differential variables.
5214 
5215    Level: beginner
5216 
5217 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
5218 @*/
5219 PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
5220 {
5221   PetscErrorCode ierr;
5222 
5223   PetscFunctionBegin;
5224   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
5225   if (vatol) {
5226     ierr = PetscObjectReference((PetscObject)vatol);CHKERRQ(ierr);
5227     ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
5228 
5229     ts->vatol = vatol;
5230   }
5231   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
5232   if (vrtol) {
5233     ierr = PetscObjectReference((PetscObject)vrtol);CHKERRQ(ierr);
5234     ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
5235 
5236     ts->vrtol = vrtol;
5237   }
5238   PetscFunctionReturn(0);
5239 }
5240 
5241 #undef __FUNCT__
5242 #define __FUNCT__ "TSGetTolerances"
5243 /*@
5244    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
5245 
5246    Logically Collective
5247 
5248    Input Arguments:
5249 .  ts - time integration context
5250 
5251    Output Arguments:
5252 +  atol - scalar absolute tolerances, NULL to ignore
5253 .  vatol - vector of absolute tolerances, NULL to ignore
5254 .  rtol - scalar relative tolerances, NULL to ignore
5255 -  vrtol - vector of relative tolerances, NULL to ignore
5256 
5257    Level: beginner
5258 
5259 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
5260 @*/
5261 PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
5262 {
5263   PetscFunctionBegin;
5264   if (atol)  *atol  = ts->atol;
5265   if (vatol) *vatol = ts->vatol;
5266   if (rtol)  *rtol  = ts->rtol;
5267   if (vrtol) *vrtol = ts->vrtol;
5268   PetscFunctionReturn(0);
5269 }
5270 
5271 #undef __FUNCT__
5272 #define __FUNCT__ "TSErrorWeightedNorm2"
5273 /*@
5274    TSErrorWeightedNorm2 - compute a weighted 2-norm of the difference between two state vectors
5275 
5276    Collective on TS
5277 
5278    Input Arguments:
5279 +  ts - time stepping context
5280 .  U - state vector, usually ts->vec_sol
5281 -  Y - state vector to be compared to U
5282 
5283    Output Arguments:
5284 .  norm - weighted norm, a value of 1.0 is considered small
5285 
5286    Level: developer
5287 
5288 .seealso: TSErrorWeightedNorm(), TSErrorWeightedNormInfinity()
5289 @*/
5290 PetscErrorCode TSErrorWeightedNorm2(TS ts,Vec U,Vec Y,PetscReal *norm)
5291 {
5292   PetscErrorCode    ierr;
5293   PetscInt          i,n,N,rstart;
5294   const PetscScalar *u,*y;
5295   PetscReal         sum,gsum;
5296   PetscReal         tol;
5297 
5298   PetscFunctionBegin;
5299   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5300   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
5301   PetscValidHeaderSpecific(Y,VEC_CLASSID,3);
5302   PetscValidType(U,2);
5303   PetscValidType(Y,3);
5304   PetscCheckSameComm(U,2,Y,3);
5305   PetscValidPointer(norm,4);
5306   if (U == Y) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_IDN,"U and Y cannot be the same vector");
5307 
5308   ierr = VecGetSize(U,&N);CHKERRQ(ierr);
5309   ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr);
5310   ierr = VecGetOwnershipRange(U,&rstart,NULL);CHKERRQ(ierr);
5311   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
5312   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
5313   sum  = 0.;
5314   if (ts->vatol && ts->vrtol) {
5315     const PetscScalar *atol,*rtol;
5316     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
5317     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
5318     for (i=0; i<n; i++) {
5319       tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5320       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
5321     }
5322     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
5323     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
5324   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5325     const PetscScalar *atol;
5326     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
5327     for (i=0; i<n; i++) {
5328       tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5329       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
5330     }
5331     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
5332   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5333     const PetscScalar *rtol;
5334     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
5335     for (i=0; i<n; i++) {
5336       tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5337       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
5338     }
5339     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
5340   } else {                      /* scalar atol, scalar rtol */
5341     for (i=0; i<n; i++) {
5342       tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5343       sum += PetscSqr(PetscAbsScalar(y[i] - u[i]) / tol);
5344     }
5345   }
5346   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
5347   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
5348 
5349   ierr  = MPIU_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
5350   *norm = PetscSqrtReal(gsum / N);
5351 
5352   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5353   PetscFunctionReturn(0);
5354 }
5355 
5356 #undef __FUNCT__
5357 #define __FUNCT__ "TSErrorWeightedNormInfinity"
5358 /*@
5359    TSErrorWeightedNormInfinity - compute a weighted infinity-norm of the difference between two state vectors
5360 
5361    Collective on TS
5362 
5363    Input Arguments:
5364 +  ts - time stepping context
5365 .  U - state vector, usually ts->vec_sol
5366 -  Y - state vector to be compared to U
5367 
5368    Output Arguments:
5369 .  norm - weighted norm, a value of 1.0 is considered small
5370 
5371    Level: developer
5372 
5373 .seealso: TSErrorWeightedNorm(), TSErrorWeightedNorm2()
5374 @*/
5375 PetscErrorCode TSErrorWeightedNormInfinity(TS ts,Vec U,Vec Y,PetscReal *norm)
5376 {
5377   PetscErrorCode    ierr;
5378   PetscInt          i,n,N,rstart,k;
5379   const PetscScalar *u,*y;
5380   PetscReal         max,gmax;
5381   PetscReal         tol;
5382 
5383   PetscFunctionBegin;
5384   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5385   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
5386   PetscValidHeaderSpecific(Y,VEC_CLASSID,3);
5387   PetscValidType(U,2);
5388   PetscValidType(Y,3);
5389   PetscCheckSameComm(U,2,Y,3);
5390   PetscValidPointer(norm,4);
5391   if (U == Y) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_IDN,"U and Y cannot be the same vector");
5392 
5393   ierr = VecGetSize(U,&N);CHKERRQ(ierr);
5394   ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr);
5395   ierr = VecGetOwnershipRange(U,&rstart,NULL);CHKERRQ(ierr);
5396   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
5397   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
5398   if (ts->vatol && ts->vrtol) {
5399     const PetscScalar *atol,*rtol;
5400     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
5401     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
5402     k = 0;
5403     tol = PetscRealPart(atol[k]) + PetscRealPart(rtol[k]) * PetscMax(PetscAbsScalar(u[k]),PetscAbsScalar(y[k]));
5404     max = PetscAbsScalar(y[k] - u[k]) / tol;
5405     for (i=1; i<n; i++) {
5406       tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5407       max = PetscMax(max,PetscAbsScalar(y[i] - u[i]) / tol);
5408     }
5409     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
5410     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
5411   } else if (ts->vatol) {       /* vector atol, scalar rtol */
5412     const PetscScalar *atol;
5413     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
5414     k = 0;
5415     tol = PetscRealPart(atol[k]) + ts->rtol * PetscMax(PetscAbsScalar(u[k]),PetscAbsScalar(y[k]));
5416     max = PetscAbsScalar(y[k] - u[k]) / tol;
5417     for (i=1; i<n; i++) {
5418       tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5419       max = PetscMax(max,PetscAbsScalar(y[i] - u[i]) / tol);
5420     }
5421     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
5422   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
5423     const PetscScalar *rtol;
5424     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
5425     k = 0;
5426     tol = ts->atol + PetscRealPart(rtol[k]) * PetscMax(PetscAbsScalar(u[k]),PetscAbsScalar(y[k]));
5427     max = PetscAbsScalar(y[k] - u[k]) / tol;
5428     for (i=1; i<n; i++) {
5429       tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5430       max = PetscMax(max,PetscAbsScalar(y[i] - u[i]) / tol);
5431     }
5432     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
5433   } else {                      /* scalar atol, scalar rtol */
5434     k = 0;
5435     tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(u[k]),PetscAbsScalar(y[k]));
5436     max = PetscAbsScalar(y[k] - u[k]) / tol;
5437     for (i=1; i<n; i++) {
5438       tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(u[i]),PetscAbsScalar(y[i]));
5439       max = PetscMax(max,PetscAbsScalar(y[i] - u[i]) / tol);
5440     }
5441   }
5442   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
5443   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
5444 
5445   ierr  = MPIU_Allreduce(&max,&gmax,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
5446   *norm = gmax;
5447 
5448   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
5449   PetscFunctionReturn(0);
5450 }
5451 
5452 #undef __FUNCT__
5453 #define __FUNCT__ "TSErrorWeightedNorm"
5454 /*@
5455    TSErrorWeightedNorm - compute a weighted norm of the difference between two state vectors
5456 
5457    Collective on TS
5458 
5459    Input Arguments:
5460 +  ts - time stepping context
5461 .  U - state vector, usually ts->vec_sol
5462 .  Y - state vector to be compared to U
5463 -  wnormtype - norm type, either NORM_2 or NORM_INFINITY
5464 
5465    Output Arguments:
5466 .  norm - weighted norm, a value of 1.0 is considered small
5467 
5468 
5469    Options Database Keys:
5470 .  -ts_adapt_wnormtype <wnormtype> - 2, INFINITY
5471 
5472    Level: developer
5473 
5474 .seealso: TSErrorWeightedNormInfinity(), TSErrorWeightedNorm2()
5475 @*/
5476 PetscErrorCode TSErrorWeightedNorm(TS ts,Vec U,Vec Y,NormType wnormtype,PetscReal *norm)
5477 {
5478   PetscErrorCode ierr;
5479 
5480   PetscFunctionBegin;
5481   if (wnormtype == NORM_2) {
5482     ierr = TSErrorWeightedNorm2(ts,U,Y,norm);CHKERRQ(ierr);
5483   } else if(wnormtype == NORM_INFINITY) {
5484     ierr = TSErrorWeightedNormInfinity(ts,U,Y,norm);CHKERRQ(ierr);
5485   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for norm type %s",NormTypes[wnormtype]);
5486   PetscFunctionReturn(0);
5487 }
5488 
5489 #undef __FUNCT__
5490 #define __FUNCT__ "TSSetCFLTimeLocal"
5491 /*@
5492    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
5493 
5494    Logically Collective on TS
5495 
5496    Input Arguments:
5497 +  ts - time stepping context
5498 -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)
5499 
5500    Note:
5501    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
5502 
5503    Level: intermediate
5504 
5505 .seealso: TSGetCFLTime(), TSADAPTCFL
5506 @*/
5507 PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
5508 {
5509   PetscFunctionBegin;
5510   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5511   ts->cfltime_local = cfltime;
5512   ts->cfltime       = -1.;
5513   PetscFunctionReturn(0);
5514 }
5515 
5516 #undef __FUNCT__
5517 #define __FUNCT__ "TSGetCFLTime"
5518 /*@
5519    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
5520 
5521    Collective on TS
5522 
5523    Input Arguments:
5524 .  ts - time stepping context
5525 
5526    Output Arguments:
5527 .  cfltime - maximum stable time step for forward Euler
5528 
5529    Level: advanced
5530 
5531 .seealso: TSSetCFLTimeLocal()
5532 @*/
5533 PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
5534 {
5535   PetscErrorCode ierr;
5536 
5537   PetscFunctionBegin;
5538   if (ts->cfltime < 0) {
5539     ierr = MPIU_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
5540   }
5541   *cfltime = ts->cfltime;
5542   PetscFunctionReturn(0);
5543 }
5544 
5545 #undef __FUNCT__
5546 #define __FUNCT__ "TSVISetVariableBounds"
5547 /*@
5548    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
5549 
5550    Input Parameters:
5551 .  ts   - the TS context.
5552 .  xl   - lower bound.
5553 .  xu   - upper bound.
5554 
5555    Notes:
5556    If this routine is not called then the lower and upper bounds are set to
5557    PETSC_NINFINITY and PETSC_INFINITY respectively during SNESSetUp().
5558 
5559    Level: advanced
5560 
5561 @*/
5562 PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
5563 {
5564   PetscErrorCode ierr;
5565   SNES           snes;
5566 
5567   PetscFunctionBegin;
5568   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
5569   ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr);
5570   PetscFunctionReturn(0);
5571 }
5572 
5573 #if defined(PETSC_HAVE_MATLAB_ENGINE)
5574 #include <mex.h>
5575 
5576 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
5577 
5578 #undef __FUNCT__
5579 #define __FUNCT__ "TSComputeFunction_Matlab"
5580 /*
5581    TSComputeFunction_Matlab - Calls the function that has been set with
5582                          TSSetFunctionMatlab().
5583 
5584    Collective on TS
5585 
5586    Input Parameters:
5587 +  snes - the TS context
5588 -  u - input vector
5589 
5590    Output Parameter:
5591 .  y - function vector, as set by TSSetFunction()
5592 
5593    Notes:
5594    TSComputeFunction() is typically used within nonlinear solvers
5595    implementations, so most users would not generally call this routine
5596    themselves.
5597 
5598    Level: developer
5599 
5600 .keywords: TS, nonlinear, compute, function
5601 
5602 .seealso: TSSetFunction(), TSGetFunction()
5603 */
5604 PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec u,Vec udot,Vec y, void *ctx)
5605 {
5606   PetscErrorCode  ierr;
5607   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
5608   int             nlhs  = 1,nrhs = 7;
5609   mxArray         *plhs[1],*prhs[7];
5610   long long int   lx = 0,lxdot = 0,ly = 0,ls = 0;
5611 
5612   PetscFunctionBegin;
5613   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
5614   PetscValidHeaderSpecific(u,VEC_CLASSID,3);
5615   PetscValidHeaderSpecific(udot,VEC_CLASSID,4);
5616   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
5617   PetscCheckSameComm(snes,1,u,3);
5618   PetscCheckSameComm(snes,1,y,5);
5619 
5620   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
5621   ierr = PetscMemcpy(&lx,&u,sizeof(u));CHKERRQ(ierr);
5622   ierr = PetscMemcpy(&lxdot,&udot,sizeof(udot));CHKERRQ(ierr);
5623   ierr = PetscMemcpy(&ly,&y,sizeof(u));CHKERRQ(ierr);
5624 
5625   prhs[0] =  mxCreateDoubleScalar((double)ls);
5626   prhs[1] =  mxCreateDoubleScalar(time);
5627   prhs[2] =  mxCreateDoubleScalar((double)lx);
5628   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
5629   prhs[4] =  mxCreateDoubleScalar((double)ly);
5630   prhs[5] =  mxCreateString(sctx->funcname);
5631   prhs[6] =  sctx->ctx;
5632   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
5633   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
5634   mxDestroyArray(prhs[0]);
5635   mxDestroyArray(prhs[1]);
5636   mxDestroyArray(prhs[2]);
5637   mxDestroyArray(prhs[3]);
5638   mxDestroyArray(prhs[4]);
5639   mxDestroyArray(prhs[5]);
5640   mxDestroyArray(plhs[0]);
5641   PetscFunctionReturn(0);
5642 }
5643 
5644 
5645 #undef __FUNCT__
5646 #define __FUNCT__ "TSSetFunctionMatlab"
5647 /*
5648    TSSetFunctionMatlab - Sets the function evaluation routine and function
5649    vector for use by the TS routines in solving ODEs
5650    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
5651 
5652    Logically Collective on TS
5653 
5654    Input Parameters:
5655 +  ts - the TS context
5656 -  func - function evaluation routine
5657 
5658    Calling sequence of func:
5659 $    func (TS ts,PetscReal time,Vec u,Vec udot,Vec f,void *ctx);
5660 
5661    Level: beginner
5662 
5663 .keywords: TS, nonlinear, set, function
5664 
5665 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
5666 */
5667 PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
5668 {
5669   PetscErrorCode  ierr;
5670   TSMatlabContext *sctx;
5671 
5672   PetscFunctionBegin;
5673   /* currently sctx is memory bleed */
5674   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
5675   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
5676   /*
5677      This should work, but it doesn't
5678   sctx->ctx = ctx;
5679   mexMakeArrayPersistent(sctx->ctx);
5680   */
5681   sctx->ctx = mxDuplicateArray(ctx);
5682 
5683   ierr = TSSetIFunction(ts,NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
5684   PetscFunctionReturn(0);
5685 }
5686 
5687 #undef __FUNCT__
5688 #define __FUNCT__ "TSComputeJacobian_Matlab"
5689 /*
5690    TSComputeJacobian_Matlab - Calls the function that has been set with
5691                          TSSetJacobianMatlab().
5692 
5693    Collective on TS
5694 
5695    Input Parameters:
5696 +  ts - the TS context
5697 .  u - input vector
5698 .  A, B - the matrices
5699 -  ctx - user context
5700 
5701    Level: developer
5702 
5703 .keywords: TS, nonlinear, compute, function
5704 
5705 .seealso: TSSetFunction(), TSGetFunction()
5706 @*/
5707 PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec u,Vec udot,PetscReal shift,Mat A,Mat B,void *ctx)
5708 {
5709   PetscErrorCode  ierr;
5710   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
5711   int             nlhs  = 2,nrhs = 9;
5712   mxArray         *plhs[2],*prhs[9];
5713   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
5714 
5715   PetscFunctionBegin;
5716   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5717   PetscValidHeaderSpecific(u,VEC_CLASSID,3);
5718 
5719   /* call Matlab function in ctx with arguments u and y */
5720 
5721   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
5722   ierr = PetscMemcpy(&lx,&u,sizeof(u));CHKERRQ(ierr);
5723   ierr = PetscMemcpy(&lxdot,&udot,sizeof(u));CHKERRQ(ierr);
5724   ierr = PetscMemcpy(&lA,A,sizeof(u));CHKERRQ(ierr);
5725   ierr = PetscMemcpy(&lB,B,sizeof(u));CHKERRQ(ierr);
5726 
5727   prhs[0] =  mxCreateDoubleScalar((double)ls);
5728   prhs[1] =  mxCreateDoubleScalar((double)time);
5729   prhs[2] =  mxCreateDoubleScalar((double)lx);
5730   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
5731   prhs[4] =  mxCreateDoubleScalar((double)shift);
5732   prhs[5] =  mxCreateDoubleScalar((double)lA);
5733   prhs[6] =  mxCreateDoubleScalar((double)lB);
5734   prhs[7] =  mxCreateString(sctx->funcname);
5735   prhs[8] =  sctx->ctx;
5736   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
5737   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
5738   mxDestroyArray(prhs[0]);
5739   mxDestroyArray(prhs[1]);
5740   mxDestroyArray(prhs[2]);
5741   mxDestroyArray(prhs[3]);
5742   mxDestroyArray(prhs[4]);
5743   mxDestroyArray(prhs[5]);
5744   mxDestroyArray(prhs[6]);
5745   mxDestroyArray(prhs[7]);
5746   mxDestroyArray(plhs[0]);
5747   mxDestroyArray(plhs[1]);
5748   PetscFunctionReturn(0);
5749 }
5750 
5751 
5752 #undef __FUNCT__
5753 #define __FUNCT__ "TSSetJacobianMatlab"
5754 /*
5755    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
5756    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
5757 
5758    Logically Collective on TS
5759 
5760    Input Parameters:
5761 +  ts - the TS context
5762 .  A,B - Jacobian matrices
5763 .  func - function evaluation routine
5764 -  ctx - user context
5765 
5766    Calling sequence of func:
5767 $    flag = func (TS ts,PetscReal time,Vec u,Vec udot,Mat A,Mat B,void *ctx);
5768 
5769 
5770    Level: developer
5771 
5772 .keywords: TS, nonlinear, set, function
5773 
5774 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
5775 */
5776 PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
5777 {
5778   PetscErrorCode  ierr;
5779   TSMatlabContext *sctx;
5780 
5781   PetscFunctionBegin;
5782   /* currently sctx is memory bleed */
5783   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
5784   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
5785   /*
5786      This should work, but it doesn't
5787   sctx->ctx = ctx;
5788   mexMakeArrayPersistent(sctx->ctx);
5789   */
5790   sctx->ctx = mxDuplicateArray(ctx);
5791 
5792   ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
5793   PetscFunctionReturn(0);
5794 }
5795 
5796 #undef __FUNCT__
5797 #define __FUNCT__ "TSMonitor_Matlab"
5798 /*
5799    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
5800 
5801    Collective on TS
5802 
5803 .seealso: TSSetFunction(), TSGetFunction()
5804 @*/
5805 PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec u, void *ctx)
5806 {
5807   PetscErrorCode  ierr;
5808   TSMatlabContext *sctx = (TSMatlabContext*)ctx;
5809   int             nlhs  = 1,nrhs = 6;
5810   mxArray         *plhs[1],*prhs[6];
5811   long long int   lx = 0,ls = 0;
5812 
5813   PetscFunctionBegin;
5814   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5815   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
5816 
5817   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
5818   ierr = PetscMemcpy(&lx,&u,sizeof(u));CHKERRQ(ierr);
5819 
5820   prhs[0] =  mxCreateDoubleScalar((double)ls);
5821   prhs[1] =  mxCreateDoubleScalar((double)it);
5822   prhs[2] =  mxCreateDoubleScalar((double)time);
5823   prhs[3] =  mxCreateDoubleScalar((double)lx);
5824   prhs[4] =  mxCreateString(sctx->funcname);
5825   prhs[5] =  sctx->ctx;
5826   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
5827   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
5828   mxDestroyArray(prhs[0]);
5829   mxDestroyArray(prhs[1]);
5830   mxDestroyArray(prhs[2]);
5831   mxDestroyArray(prhs[3]);
5832   mxDestroyArray(prhs[4]);
5833   mxDestroyArray(plhs[0]);
5834   PetscFunctionReturn(0);
5835 }
5836 
5837 
5838 #undef __FUNCT__
5839 #define __FUNCT__ "TSMonitorSetMatlab"
5840 /*
5841    TSMonitorSetMatlab - Sets the monitor function from Matlab
5842 
5843    Level: developer
5844 
5845 .keywords: TS, nonlinear, set, function
5846 
5847 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
5848 */
5849 PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
5850 {
5851   PetscErrorCode  ierr;
5852   TSMatlabContext *sctx;
5853 
5854   PetscFunctionBegin;
5855   /* currently sctx is memory bleed */
5856   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
5857   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
5858   /*
5859      This should work, but it doesn't
5860   sctx->ctx = ctx;
5861   mexMakeArrayPersistent(sctx->ctx);
5862   */
5863   sctx->ctx = mxDuplicateArray(ctx);
5864 
5865   ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,NULL);CHKERRQ(ierr);
5866   PetscFunctionReturn(0);
5867 }
5868 #endif
5869 
5870 #undef __FUNCT__
5871 #define __FUNCT__ "TSMonitorLGSolution"
5872 /*@C
5873    TSMonitorLGSolution - Monitors progress of the TS solvers by plotting each component of the solution vector
5874        in a time based line graph
5875 
5876    Collective on TS
5877 
5878    Input Parameters:
5879 +  ts - the TS context
5880 .  step - current time-step
5881 .  ptime - current time
5882 .  u - current solution
5883 -  dctx - the TSMonitorLGCtx object that contains all the options for the monitoring, this is created with TSMonitorLGCtxCreate()
5884 
5885    Options Database:
5886 .   -ts_monitor_lg_solution_variables
5887 
5888    Level: intermediate
5889 
5890     Notes: each process in a parallel run displays its component solutions in a separate window
5891 
5892 .keywords: TS,  vector, monitor, view
5893 
5894 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGCtxCreate(), TSMonitorLGCtxSetVariableNames(), TSMonitorLGCtxGetVariableNames(),
5895            TSMonitorLGSetVariableNames(), TSMonitorLGGetVariableNames(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetDisplayVariables(),
5896            TSMonitorLGCtxSetTransform(), TSMonitorLGSetTransform(), TSMonitorLGError(), TSMonitorLGSNESIterations(), TSMonitorLGKSPIterations(),
5897            TSMonitorEnvelopeCtxCreate(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxDestroy(), TSMonitorEnvelop()
5898 @*/
5899 PetscErrorCode  TSMonitorLGSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx)
5900 {
5901   PetscErrorCode    ierr;
5902   TSMonitorLGCtx    ctx = (TSMonitorLGCtx)dctx;
5903   const PetscScalar *yy;
5904   PetscInt          dim;
5905   Vec               v;
5906 
5907   PetscFunctionBegin;
5908   if (!step) {
5909     PetscDrawAxis axis;
5910     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
5911     ierr = PetscDrawAxisSetLabels(axis,"Solution as function of time","Time","Solution");CHKERRQ(ierr);
5912     if (ctx->names && !ctx->displaynames) {
5913       char      **displaynames;
5914       PetscBool flg;
5915 
5916       ierr = VecGetLocalSize(u,&dim);CHKERRQ(ierr);
5917       ierr = PetscMalloc((dim+1)*sizeof(char*),&displaynames);CHKERRQ(ierr);
5918       ierr = PetscMemzero(displaynames,(dim+1)*sizeof(char*));CHKERRQ(ierr);
5919       ierr = PetscOptionsGetStringArray(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_monitor_lg_solution_variables",displaynames,&dim,&flg);CHKERRQ(ierr);
5920       if (flg) {
5921         ierr = TSMonitorLGCtxSetDisplayVariables(ctx,(const char *const *)displaynames);CHKERRQ(ierr);
5922       }
5923       ierr = PetscStrArrayDestroy(&displaynames);CHKERRQ(ierr);
5924     }
5925     if (ctx->displaynames) {
5926       ierr = PetscDrawLGSetDimension(ctx->lg,ctx->ndisplayvariables);CHKERRQ(ierr);
5927       ierr = PetscDrawLGSetLegend(ctx->lg,(const char *const *)ctx->displaynames);CHKERRQ(ierr);
5928     } else if (ctx->names) {
5929       ierr = VecGetLocalSize(u,&dim);CHKERRQ(ierr);
5930       ierr = PetscDrawLGSetDimension(ctx->lg,dim);CHKERRQ(ierr);
5931       ierr = PetscDrawLGSetLegend(ctx->lg,(const char *const *)ctx->names);CHKERRQ(ierr);
5932     } else {
5933       ierr = VecGetLocalSize(u,&dim);CHKERRQ(ierr);
5934       ierr = PetscDrawLGSetDimension(ctx->lg,dim);CHKERRQ(ierr);
5935     }
5936     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
5937   }
5938   if (ctx->transform) {
5939     ierr = (*ctx->transform)(ctx->transformctx,u,&v);CHKERRQ(ierr);
5940   } else {
5941     v = u;
5942   }
5943   ierr = VecGetArrayRead(v,&yy);CHKERRQ(ierr);
5944 #if defined(PETSC_USE_COMPLEX)
5945   {
5946     PetscReal *yreal;
5947     PetscInt  i,n;
5948     ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
5949     ierr = PetscMalloc1(n,&yreal);CHKERRQ(ierr);
5950     for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
5951     ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);CHKERRQ(ierr);
5952     ierr = PetscFree(yreal);CHKERRQ(ierr);
5953   }
5954 #else
5955   if (ctx->displaynames) {
5956     PetscInt i;
5957     for (i=0; i<ctx->ndisplayvariables; i++) {
5958       ctx->displayvalues[i] = yy[ctx->displayvariables[i]];
5959     }
5960     ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,ctx->displayvalues);CHKERRQ(ierr);
5961   } else {
5962     ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);CHKERRQ(ierr);
5963   }
5964 #endif
5965   ierr = VecRestoreArrayRead(v,&yy);CHKERRQ(ierr);
5966   if (ctx->transform) {
5967     ierr = VecDestroy(&v);CHKERRQ(ierr);
5968   }
5969   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
5970     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
5971   }
5972   PetscFunctionReturn(0);
5973 }
5974 
5975 
5976 #undef __FUNCT__
5977 #define __FUNCT__ "TSMonitorLGSetVariableNames"
5978 /*@C
5979    TSMonitorLGSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
5980 
5981    Collective on TS
5982 
5983    Input Parameters:
5984 +  ts - the TS context
5985 -  names - the names of the components, final string must be NULL
5986 
5987    Level: intermediate
5988 
5989    Notes: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored
5990 
5991 .keywords: TS,  vector, monitor, view
5992 
5993 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables(), TSMonitorLGCtxSetVariableNames()
5994 @*/
5995 PetscErrorCode  TSMonitorLGSetVariableNames(TS ts,const char * const *names)
5996 {
5997   PetscErrorCode    ierr;
5998   PetscInt          i;
5999 
6000   PetscFunctionBegin;
6001   for (i=0; i<ts->numbermonitors; i++) {
6002     if (ts->monitor[i] == TSMonitorLGSolution) {
6003       ierr = TSMonitorLGCtxSetVariableNames((TSMonitorLGCtx)ts->monitorcontext[i],names);CHKERRQ(ierr);
6004       break;
6005     }
6006   }
6007   PetscFunctionReturn(0);
6008 }
6009 
6010 #undef __FUNCT__
6011 #define __FUNCT__ "TSMonitorLGCtxSetVariableNames"
6012 /*@C
6013    TSMonitorLGCtxSetVariableNames - Sets the name of each component in the solution vector so that it may be displayed in the plot
6014 
6015    Collective on TS
6016 
6017    Input Parameters:
6018 +  ts - the TS context
6019 -  names - the names of the components, final string must be NULL
6020 
6021    Level: intermediate
6022 
6023 .keywords: TS,  vector, monitor, view
6024 
6025 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables(), TSMonitorLGSetVariableNames()
6026 @*/
6027 PetscErrorCode  TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx ctx,const char * const *names)
6028 {
6029   PetscErrorCode    ierr;
6030 
6031   PetscFunctionBegin;
6032   ierr = PetscStrArrayDestroy(&ctx->names);CHKERRQ(ierr);
6033   ierr = PetscStrArrayallocpy(names,&ctx->names);CHKERRQ(ierr);
6034   PetscFunctionReturn(0);
6035 }
6036 
6037 #undef __FUNCT__
6038 #define __FUNCT__ "TSMonitorLGGetVariableNames"
6039 /*@C
6040    TSMonitorLGGetVariableNames - Gets the name of each component in the solution vector so that it may be displayed in the plot
6041 
6042    Collective on TS
6043 
6044    Input Parameter:
6045 .  ts - the TS context
6046 
6047    Output Parameter:
6048 .  names - the names of the components, final string must be NULL
6049 
6050    Level: intermediate
6051 
6052    Notes: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored
6053 
6054 .keywords: TS,  vector, monitor, view
6055 
6056 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables()
6057 @*/
6058 PetscErrorCode  TSMonitorLGGetVariableNames(TS ts,const char *const **names)
6059 {
6060   PetscInt       i;
6061 
6062   PetscFunctionBegin;
6063   *names = NULL;
6064   for (i=0; i<ts->numbermonitors; i++) {
6065     if (ts->monitor[i] == TSMonitorLGSolution) {
6066       TSMonitorLGCtx  ctx = (TSMonitorLGCtx) ts->monitorcontext[i];
6067       *names = (const char *const *)ctx->names;
6068       break;
6069     }
6070   }
6071   PetscFunctionReturn(0);
6072 }
6073 
6074 #undef __FUNCT__
6075 #define __FUNCT__ "TSMonitorLGCtxSetDisplayVariables"
6076 /*@C
6077    TSMonitorLGCtxSetDisplayVariables - Sets the variables that are to be display in the monitor
6078 
6079    Collective on TS
6080 
6081    Input Parameters:
6082 +  ctx - the TSMonitorLG context
6083 .  displaynames - the names of the components, final string must be NULL
6084 
6085    Level: intermediate
6086 
6087 .keywords: TS,  vector, monitor, view
6088 
6089 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames()
6090 @*/
6091 PetscErrorCode  TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx ctx,const char * const *displaynames)
6092 {
6093   PetscInt          j = 0,k;
6094   PetscErrorCode    ierr;
6095 
6096   PetscFunctionBegin;
6097   if (!ctx->names) PetscFunctionReturn(0);
6098   ierr = PetscStrArrayDestroy(&ctx->displaynames);CHKERRQ(ierr);
6099   ierr = PetscStrArrayallocpy(displaynames,&ctx->displaynames);CHKERRQ(ierr);
6100   while (displaynames[j]) j++;
6101   ctx->ndisplayvariables = j;
6102   ierr = PetscMalloc1(ctx->ndisplayvariables,&ctx->displayvariables);CHKERRQ(ierr);
6103   ierr = PetscMalloc1(ctx->ndisplayvariables,&ctx->displayvalues);CHKERRQ(ierr);
6104   j = 0;
6105   while (displaynames[j]) {
6106     k = 0;
6107     while (ctx->names[k]) {
6108       PetscBool flg;
6109       ierr = PetscStrcmp(displaynames[j],ctx->names[k],&flg);CHKERRQ(ierr);
6110       if (flg) {
6111         ctx->displayvariables[j] = k;
6112         break;
6113       }
6114       k++;
6115     }
6116     j++;
6117   }
6118   PetscFunctionReturn(0);
6119 }
6120 
6121 
6122 #undef __FUNCT__
6123 #define __FUNCT__ "TSMonitorLGSetDisplayVariables"
6124 /*@C
6125    TSMonitorLGSetDisplayVariables - Sets the variables that are to be display in the monitor
6126 
6127    Collective on TS
6128 
6129    Input Parameters:
6130 +  ts - the TS context
6131 .  displaynames - the names of the components, final string must be NULL
6132 
6133    Notes: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored
6134 
6135    Level: intermediate
6136 
6137 .keywords: TS,  vector, monitor, view
6138 
6139 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames()
6140 @*/
6141 PetscErrorCode  TSMonitorLGSetDisplayVariables(TS ts,const char * const *displaynames)
6142 {
6143   PetscInt          i;
6144   PetscErrorCode    ierr;
6145 
6146   PetscFunctionBegin;
6147   for (i=0; i<ts->numbermonitors; i++) {
6148     if (ts->monitor[i] == TSMonitorLGSolution) {
6149       ierr = TSMonitorLGCtxSetDisplayVariables((TSMonitorLGCtx)ts->monitorcontext[i],displaynames);CHKERRQ(ierr);
6150       break;
6151     }
6152   }
6153   PetscFunctionReturn(0);
6154 }
6155 
6156 #undef __FUNCT__
6157 #define __FUNCT__ "TSMonitorLGSetTransform"
6158 /*@C
6159    TSMonitorLGSetTransform - Solution vector will be transformed by provided function before being displayed
6160 
6161    Collective on TS
6162 
6163    Input Parameters:
6164 +  ts - the TS context
6165 .  transform - the transform function
6166 .  destroy - function to destroy the optional context
6167 -  ctx - optional context used by transform function
6168 
6169    Notes: If the TS object does not have a TSMonitorLGCtx associated with it then this function is ignored
6170 
6171    Level: intermediate
6172 
6173 .keywords: TS,  vector, monitor, view
6174 
6175 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames(), TSMonitorLGCtxSetTransform()
6176 @*/
6177 PetscErrorCode  TSMonitorLGSetTransform(TS ts,PetscErrorCode (*transform)(void*,Vec,Vec*),PetscErrorCode (*destroy)(void*),void *tctx)
6178 {
6179   PetscInt          i;
6180   PetscErrorCode    ierr;
6181 
6182   PetscFunctionBegin;
6183   for (i=0; i<ts->numbermonitors; i++) {
6184     if (ts->monitor[i] == TSMonitorLGSolution) {
6185       ierr = TSMonitorLGCtxSetTransform((TSMonitorLGCtx)ts->monitorcontext[i],transform,destroy,tctx);CHKERRQ(ierr);
6186     }
6187   }
6188   PetscFunctionReturn(0);
6189 }
6190 
6191 #undef __FUNCT__
6192 #define __FUNCT__ "TSMonitorLGCtxSetTransform"
6193 /*@C
6194    TSMonitorLGCtxSetTransform - Solution vector will be transformed by provided function before being displayed
6195 
6196    Collective on TSLGCtx
6197 
6198    Input Parameters:
6199 +  ts - the TS context
6200 .  transform - the transform function
6201 .  destroy - function to destroy the optional context
6202 -  ctx - optional context used by transform function
6203 
6204    Level: intermediate
6205 
6206 .keywords: TS,  vector, monitor, view
6207 
6208 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetVariableNames(), TSMonitorLGSetTransform()
6209 @*/
6210 PetscErrorCode  TSMonitorLGCtxSetTransform(TSMonitorLGCtx ctx,PetscErrorCode (*transform)(void*,Vec,Vec*),PetscErrorCode (*destroy)(void*),void *tctx)
6211 {
6212   PetscFunctionBegin;
6213   ctx->transform    = transform;
6214   ctx->transformdestroy = destroy;
6215   ctx->transformctx = tctx;
6216   PetscFunctionReturn(0);
6217 }
6218 
6219 #undef __FUNCT__
6220 #define __FUNCT__ "TSMonitorLGError"
6221 /*@C
6222    TSMonitorLGError - Monitors progress of the TS solvers by plotting each component of the solution vector
6223        in a time based line graph
6224 
6225    Collective on TS
6226 
6227    Input Parameters:
6228 +  ts - the TS context
6229 .  step - current time-step
6230 .  ptime - current time
6231 .  u - current solution
6232 -  dctx - TSMonitorLGCtx object created with TSMonitorLGCtxCreate()
6233 
6234    Level: intermediate
6235 
6236    Notes:
6237    Only for sequential solves.
6238 
6239    The user must provide the solution using TSSetSolutionFunction() to use this monitor.
6240 
6241    Options Database Keys:
6242 .  -ts_monitor_lg_error - create a graphical monitor of error history
6243 
6244 .keywords: TS,  vector, monitor, view
6245 
6246 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSSetSolutionFunction()
6247 @*/
6248 PetscErrorCode  TSMonitorLGError(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dummy)
6249 {
6250   PetscErrorCode    ierr;
6251   TSMonitorLGCtx    ctx = (TSMonitorLGCtx)dummy;
6252   const PetscScalar *yy;
6253   Vec               y;
6254   PetscInt          dim;
6255 
6256   PetscFunctionBegin;
6257   if (!step) {
6258     PetscDrawAxis axis;
6259     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
6260     ierr = PetscDrawAxisSetLabels(axis,"Error in solution as function of time","Time","Solution");CHKERRQ(ierr);
6261     ierr = VecGetLocalSize(u,&dim);CHKERRQ(ierr);
6262     ierr = PetscDrawLGSetDimension(ctx->lg,dim);CHKERRQ(ierr);
6263     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
6264   }
6265   ierr = VecDuplicate(u,&y);CHKERRQ(ierr);
6266   ierr = TSComputeSolutionFunction(ts,ptime,y);CHKERRQ(ierr);
6267   ierr = VecAXPY(y,-1.0,u);CHKERRQ(ierr);
6268   ierr = VecGetArrayRead(y,&yy);CHKERRQ(ierr);
6269 #if defined(PETSC_USE_COMPLEX)
6270   {
6271     PetscReal *yreal;
6272     PetscInt  i,n;
6273     ierr = VecGetLocalSize(y,&n);CHKERRQ(ierr);
6274     ierr = PetscMalloc1(n,&yreal);CHKERRQ(ierr);
6275     for (i=0; i<n; i++) yreal[i] = PetscRealPart(yy[i]);
6276     ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,yreal);CHKERRQ(ierr);
6277     ierr = PetscFree(yreal);CHKERRQ(ierr);
6278   }
6279 #else
6280   ierr = PetscDrawLGAddCommonPoint(ctx->lg,ptime,yy);CHKERRQ(ierr);
6281 #endif
6282   ierr = VecRestoreArrayRead(y,&yy);CHKERRQ(ierr);
6283   ierr = VecDestroy(&y);CHKERRQ(ierr);
6284   if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && ts->reason)) {
6285     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
6286   }
6287   PetscFunctionReturn(0);
6288 }
6289 
6290 #undef __FUNCT__
6291 #define __FUNCT__ "TSMonitorLGSNESIterations"
6292 PetscErrorCode TSMonitorLGSNESIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
6293 {
6294   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
6295   PetscReal      x   = ptime,y;
6296   PetscErrorCode ierr;
6297   PetscInt       its;
6298 
6299   PetscFunctionBegin;
6300   if (!n) {
6301     PetscDrawAxis axis;
6302 
6303     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
6304     ierr = PetscDrawAxisSetLabels(axis,"Nonlinear iterations as function of time","Time","SNES Iterations");CHKERRQ(ierr);
6305     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
6306 
6307     ctx->snes_its = 0;
6308   }
6309   ierr = TSGetSNESIterations(ts,&its);CHKERRQ(ierr);
6310   y    = its - ctx->snes_its;
6311   ierr = PetscDrawLGAddPoint(ctx->lg,&x,&y);CHKERRQ(ierr);
6312   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
6313     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
6314   }
6315   ctx->snes_its = its;
6316   PetscFunctionReturn(0);
6317 }
6318 
6319 #undef __FUNCT__
6320 #define __FUNCT__ "TSMonitorLGKSPIterations"
6321 PetscErrorCode TSMonitorLGKSPIterations(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
6322 {
6323   TSMonitorLGCtx ctx = (TSMonitorLGCtx) monctx;
6324   PetscReal      x   = ptime,y;
6325   PetscErrorCode ierr;
6326   PetscInt       its;
6327 
6328   PetscFunctionBegin;
6329   if (!n) {
6330     PetscDrawAxis axis;
6331 
6332     ierr = PetscDrawLGGetAxis(ctx->lg,&axis);CHKERRQ(ierr);
6333     ierr = PetscDrawAxisSetLabels(axis,"Linear iterations as function of time","Time","KSP Iterations");CHKERRQ(ierr);
6334     ierr = PetscDrawLGReset(ctx->lg);CHKERRQ(ierr);
6335 
6336     ctx->ksp_its = 0;
6337   }
6338   ierr = TSGetKSPIterations(ts,&its);CHKERRQ(ierr);
6339   y    = its - ctx->ksp_its;
6340   ierr = PetscDrawLGAddPoint(ctx->lg,&x,&y);CHKERRQ(ierr);
6341   if (((ctx->howoften > 0) && (!(n % ctx->howoften)) && (n > -1)) || ((ctx->howoften == -1) && (n == -1))) {
6342     ierr = PetscDrawLGDraw(ctx->lg);CHKERRQ(ierr);
6343   }
6344   ctx->ksp_its = its;
6345   PetscFunctionReturn(0);
6346 }
6347 
6348 #undef __FUNCT__
6349 #define __FUNCT__ "TSComputeLinearStability"
6350 /*@
6351    TSComputeLinearStability - computes the linear stability function at a point
6352 
6353    Collective on TS and Vec
6354 
6355    Input Parameters:
6356 +  ts - the TS context
6357 -  xr,xi - real and imaginary part of input arguments
6358 
6359    Output Parameters:
6360 .  yr,yi - real and imaginary part of function value
6361 
6362    Level: developer
6363 
6364 .keywords: TS, compute
6365 
6366 .seealso: TSSetRHSFunction(), TSComputeIFunction()
6367 @*/
6368 PetscErrorCode TSComputeLinearStability(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
6369 {
6370   PetscErrorCode ierr;
6371 
6372   PetscFunctionBegin;
6373   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6374   if (!ts->ops->linearstability) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Linearized stability function not provided for this method");
6375   ierr = (*ts->ops->linearstability)(ts,xr,xi,yr,yi);CHKERRQ(ierr);
6376   PetscFunctionReturn(0);
6377 }
6378 
6379 /* ------------------------------------------------------------------------*/
6380 #undef __FUNCT__
6381 #define __FUNCT__ "TSMonitorEnvelopeCtxCreate"
6382 /*@C
6383    TSMonitorEnvelopeCtxCreate - Creates a context for use with TSMonitorEnvelope()
6384 
6385    Collective on TS
6386 
6387    Input Parameters:
6388 .  ts  - the ODE solver object
6389 
6390    Output Parameter:
6391 .  ctx - the context
6392 
6393    Level: intermediate
6394 
6395 .keywords: TS, monitor, line graph, residual, seealso
6396 
6397 .seealso: TSMonitorLGTimeStep(), TSMonitorSet(), TSMonitorLGSolution(), TSMonitorLGError()
6398 
6399 @*/
6400 PetscErrorCode  TSMonitorEnvelopeCtxCreate(TS ts,TSMonitorEnvelopeCtx *ctx)
6401 {
6402   PetscErrorCode ierr;
6403 
6404   PetscFunctionBegin;
6405   ierr = PetscNew(ctx);CHKERRQ(ierr);
6406   PetscFunctionReturn(0);
6407 }
6408 
6409 #undef __FUNCT__
6410 #define __FUNCT__ "TSMonitorEnvelope"
6411 /*@C
6412    TSMonitorEnvelope - Monitors the maximum and minimum value of each component of the solution
6413 
6414    Collective on TS
6415 
6416    Input Parameters:
6417 +  ts - the TS context
6418 .  step - current time-step
6419 .  ptime - current time
6420 .  u  - current solution
6421 -  dctx - the envelope context
6422 
6423    Options Database:
6424 .  -ts_monitor_envelope
6425 
6426    Level: intermediate
6427 
6428    Notes: after a solve you can use TSMonitorEnvelopeGetBounds() to access the envelope
6429 
6430 .keywords: TS,  vector, monitor, view
6431 
6432 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorEnvelopeGetBounds(), TSMonitorEnvelopeCtxCreate()
6433 @*/
6434 PetscErrorCode  TSMonitorEnvelope(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx)
6435 {
6436   PetscErrorCode       ierr;
6437   TSMonitorEnvelopeCtx ctx = (TSMonitorEnvelopeCtx)dctx;
6438 
6439   PetscFunctionBegin;
6440   if (!ctx->max) {
6441     ierr = VecDuplicate(u,&ctx->max);CHKERRQ(ierr);
6442     ierr = VecDuplicate(u,&ctx->min);CHKERRQ(ierr);
6443     ierr = VecCopy(u,ctx->max);CHKERRQ(ierr);
6444     ierr = VecCopy(u,ctx->min);CHKERRQ(ierr);
6445   } else {
6446     ierr = VecPointwiseMax(ctx->max,u,ctx->max);CHKERRQ(ierr);
6447     ierr = VecPointwiseMin(ctx->min,u,ctx->min);CHKERRQ(ierr);
6448   }
6449   PetscFunctionReturn(0);
6450 }
6451 
6452 
6453 #undef __FUNCT__
6454 #define __FUNCT__ "TSMonitorEnvelopeGetBounds"
6455 /*@C
6456    TSMonitorEnvelopeGetBounds - Gets the bounds for the components of the solution
6457 
6458    Collective on TS
6459 
6460    Input Parameter:
6461 .  ts - the TS context
6462 
6463    Output Parameter:
6464 +  max - the maximum values
6465 -  min - the minimum values
6466 
6467    Notes: If the TS does not have a TSMonitorEnvelopeCtx associated with it then this function is ignored
6468 
6469    Level: intermediate
6470 
6471 .keywords: TS,  vector, monitor, view
6472 
6473 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorLGSetDisplayVariables()
6474 @*/
6475 PetscErrorCode  TSMonitorEnvelopeGetBounds(TS ts,Vec *max,Vec *min)
6476 {
6477   PetscInt i;
6478 
6479   PetscFunctionBegin;
6480   if (max) *max = NULL;
6481   if (min) *min = NULL;
6482   for (i=0; i<ts->numbermonitors; i++) {
6483     if (ts->monitor[i] == TSMonitorEnvelope) {
6484       TSMonitorEnvelopeCtx  ctx = (TSMonitorEnvelopeCtx) ts->monitorcontext[i];
6485       if (max) *max = ctx->max;
6486       if (min) *min = ctx->min;
6487       break;
6488     }
6489   }
6490   PetscFunctionReturn(0);
6491 }
6492 
6493 #undef __FUNCT__
6494 #define __FUNCT__ "TSMonitorEnvelopeCtxDestroy"
6495 /*@C
6496    TSMonitorEnvelopeCtxDestroy - Destroys a context that was created  with TSMonitorEnvelopeCtxCreate().
6497 
6498    Collective on TSMonitorEnvelopeCtx
6499 
6500    Input Parameter:
6501 .  ctx - the monitor context
6502 
6503    Level: intermediate
6504 
6505 .keywords: TS, monitor, line graph, destroy
6506 
6507 .seealso: TSMonitorLGCtxCreate(),  TSMonitorSet(), TSMonitorLGTimeStep()
6508 @*/
6509 PetscErrorCode  TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *ctx)
6510 {
6511   PetscErrorCode ierr;
6512 
6513   PetscFunctionBegin;
6514   ierr = VecDestroy(&(*ctx)->min);CHKERRQ(ierr);
6515   ierr = VecDestroy(&(*ctx)->max);CHKERRQ(ierr);
6516   ierr = PetscFree(*ctx);CHKERRQ(ierr);
6517   PetscFunctionReturn(0);
6518 }
6519 
6520 #undef __FUNCT__
6521 #define __FUNCT__ "TSRollBack"
6522 /*@
6523    TSRollBack - Rolls back one time step
6524 
6525    Collective on TS
6526 
6527    Input Parameter:
6528 .  ts - the TS context obtained from TSCreate()
6529 
6530    Level: advanced
6531 
6532 .keywords: TS, timestep, rollback
6533 
6534 .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage(), TSInterpolate()
6535 @*/
6536 PetscErrorCode  TSRollBack(TS ts)
6537 {
6538   PetscErrorCode ierr;
6539 
6540   PetscFunctionBegin;
6541   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
6542 
6543   if (!ts->ops->rollback) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRollBack not implemented for type '%s'",((PetscObject)ts)->type_name);
6544   ierr = (*ts->ops->rollback)(ts);CHKERRQ(ierr);
6545   ts->time_step = ts->ptime - ts->ptime_prev;
6546   ts->ptime = ts->ptime_prev;
6547   ts->steprollback = PETSC_TRUE; /* Flag to indicate that the step is rollbacked */
6548   PetscFunctionReturn(0);
6549 }
6550 
6551 #undef __FUNCT__
6552 #define __FUNCT__ "TSGetStages"
6553 /*@
6554    TSGetStages - Get the number of stages and stage values
6555 
6556    Input Parameter:
6557 .  ts - the TS context obtained from TSCreate()
6558 
6559    Level: advanced
6560 
6561 .keywords: TS, getstages
6562 
6563 .seealso: TSCreate()
6564 @*/
6565 PetscErrorCode  TSGetStages(TS ts,PetscInt *ns, Vec **Y)
6566 {
6567   PetscErrorCode ierr;
6568 
6569   PetscFunctionBegin;
6570   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
6571   PetscValidPointer(ns,2);
6572 
6573   if (!ts->ops->getstages) *ns=0;
6574   else {
6575     ierr = (*ts->ops->getstages)(ts,ns,Y);CHKERRQ(ierr);
6576   }
6577   PetscFunctionReturn(0);
6578 }
6579 
6580 #undef __FUNCT__
6581 #define __FUNCT__ "TSComputeIJacobianDefaultColor"
6582 /*@C
6583   TSComputeIJacobianDefaultColor - Computes the Jacobian using finite differences and coloring to exploit matrix sparsity.
6584 
6585   Collective on SNES
6586 
6587   Input Parameters:
6588 + ts - the TS context
6589 . t - current timestep
6590 . U - state vector
6591 . Udot - time derivative of state vector
6592 . shift - shift to apply, see note below
6593 - ctx - an optional user context
6594 
6595   Output Parameters:
6596 + J - Jacobian matrix (not altered in this routine)
6597 - B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
6598 
6599   Level: intermediate
6600 
6601   Notes:
6602   If F(t,U,Udot)=0 is the DAE, the required Jacobian is
6603 
6604   dF/dU + shift*dF/dUdot
6605 
6606   Most users should not need to explicitly call this routine, as it
6607   is used internally within the nonlinear solvers.
6608 
6609   This will first try to get the coloring from the DM.  If the DM type has no coloring
6610   routine, then it will try to get the coloring from the matrix.  This requires that the
6611   matrix have nonzero entries precomputed.
6612 
6613 .keywords: TS, finite differences, Jacobian, coloring, sparse
6614 .seealso: TSSetIJacobian(), MatFDColoringCreate(), MatFDColoringSetFunction()
6615 @*/
6616 PetscErrorCode TSComputeIJacobianDefaultColor(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat J,Mat B,void *ctx)
6617 {
6618   SNES           snes;
6619   MatFDColoring  color;
6620   PetscBool      hascolor, matcolor = PETSC_FALSE;
6621   PetscErrorCode ierr;
6622 
6623   PetscFunctionBegin;
6624   ierr = PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject) ts)->prefix, "-ts_fd_color_use_mat", &matcolor, NULL);CHKERRQ(ierr);
6625   ierr = PetscObjectQuery((PetscObject) B, "TSMatFDColoring", (PetscObject *) &color);CHKERRQ(ierr);
6626   if (!color) {
6627     DM         dm;
6628     ISColoring iscoloring;
6629 
6630     ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
6631     ierr = DMHasColoring(dm, &hascolor);CHKERRQ(ierr);
6632     if (hascolor && !matcolor) {
6633       ierr = DMCreateColoring(dm, IS_COLORING_GLOBAL, &iscoloring);CHKERRQ(ierr);
6634       ierr = MatFDColoringCreate(B, iscoloring, &color);CHKERRQ(ierr);
6635       ierr = MatFDColoringSetFunction(color, (PetscErrorCode (*)(void)) SNESTSFormFunction, (void *) ts);CHKERRQ(ierr);
6636       ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
6637       ierr = MatFDColoringSetUp(B, iscoloring, color);CHKERRQ(ierr);
6638       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
6639     } else {
6640       MatColoring mc;
6641 
6642       ierr = MatColoringCreate(B, &mc);CHKERRQ(ierr);
6643       ierr = MatColoringSetDistance(mc, 2);CHKERRQ(ierr);
6644       ierr = MatColoringSetType(mc, MATCOLORINGSL);CHKERRQ(ierr);
6645       ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
6646       ierr = MatColoringApply(mc, &iscoloring);CHKERRQ(ierr);
6647       ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
6648       ierr = MatFDColoringCreate(B, iscoloring, &color);CHKERRQ(ierr);
6649       ierr = MatFDColoringSetFunction(color, (PetscErrorCode (*)(void)) SNESTSFormFunction, (void *) ts);CHKERRQ(ierr);
6650       ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
6651       ierr = MatFDColoringSetUp(B, iscoloring, color);CHKERRQ(ierr);
6652       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
6653     }
6654     ierr = PetscObjectCompose((PetscObject) B, "TSMatFDColoring", (PetscObject) color);CHKERRQ(ierr);
6655     ierr = PetscObjectDereference((PetscObject) color);CHKERRQ(ierr);
6656   }
6657   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
6658   ierr = MatFDColoringApply(B, color, U, snes);CHKERRQ(ierr);
6659   if (J != B) {
6660     ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6661     ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6662   }
6663   PetscFunctionReturn(0);
6664 }
6665 
6666 #undef __FUNCT__
6667 #define __FUNCT__ "TSSetFunctionDomainError"
6668 /*@
6669     TSSetFunctionDomainError - Set the function testing if the current state vector is valid
6670 
6671     Input Parameters:
6672     ts - the TS context
6673     func - function called within TSFunctionDomainError
6674 
6675     Level: intermediate
6676 
6677 .keywords: TS, state, domain
6678 .seealso: TSAdaptCheckStage(), TSFunctionDomainError()
6679 @*/
6680 
6681 PetscErrorCode TSSetFunctionDomainError(TS ts, PetscErrorCode (*func)(TS,PetscReal,Vec,PetscBool*))
6682 {
6683   PetscFunctionBegin;
6684   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
6685   ts->functiondomainerror = func;
6686   PetscFunctionReturn(0);
6687 }
6688 
6689 #undef __FUNCT__
6690 #define __FUNCT__ "TSFunctionDomainError"
6691 /*@
6692     TSFunctionDomainError - Check if the current state is valid
6693 
6694     Input Parameters:
6695     ts - the TS context
6696     stagetime - time of the simulation
6697     Y - state vector to check.
6698 
6699     Output Parameter:
6700     accept - Set to PETSC_FALSE if the current state vector is valid.
6701 
6702     Note:
6703     This function should be used to ensure the state is in a valid part of the space.
6704     For example, one can ensure here all values are positive.
6705 @*/
6706 PetscErrorCode TSFunctionDomainError(TS ts,PetscReal stagetime,Vec Y,PetscBool* accept)
6707 {
6708   PetscErrorCode ierr;
6709 
6710   PetscFunctionBegin;
6711 
6712   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6713   *accept = PETSC_TRUE;
6714   if (ts->functiondomainerror) {
6715     PetscStackCallStandard((*ts->functiondomainerror),(ts,stagetime,Y,accept));
6716   }
6717   PetscFunctionReturn(0);
6718 }
6719 
6720 #undef  __FUNCT__
6721 #define __FUNCT__ "TSClone"
6722 /*@C
6723   TSClone - This function clones a time step object.
6724 
6725   Collective on MPI_Comm
6726 
6727   Input Parameter:
6728 . tsin    - The input TS
6729 
6730   Output Parameter:
6731 . tsout   - The output TS (cloned)
6732 
6733   Notes:
6734   This function is used to create a clone of a TS object. It is used in ARKIMEX for initializing the slope for first stage explicit methods. It will likely be replaced in the future with a mechanism of switching methods on the fly.
6735 
6736   When using TSDestroy() on a clone the user has to first reset the correct TS reference in the embedded SNES object: e.g.: by running SNES snes_dup=NULL; TSGetSNES(ts,&snes_dup); ierr = TSSetSNES(ts,snes_dup);
6737 
6738   Level: developer
6739 
6740 .keywords: TS, clone
6741 .seealso: TSCreate(), TSSetType(), TSSetUp(), TSDestroy(), TSSetProblemType()
6742 @*/
6743 PetscErrorCode  TSClone(TS tsin, TS *tsout)
6744 {
6745   TS             t;
6746   PetscErrorCode ierr;
6747   SNES           snes_start;
6748   DM             dm;
6749   TSType         type;
6750 
6751   PetscFunctionBegin;
6752   PetscValidPointer(tsin,1);
6753   *tsout = NULL;
6754 
6755   ierr = PetscHeaderCreate(t, TS_CLASSID, "TS", "Time stepping", "TS", PetscObjectComm((PetscObject)tsin), TSDestroy, TSView);CHKERRQ(ierr);
6756 
6757   /* General TS description */
6758   t->numbermonitors    = 0;
6759   t->setupcalled       = 0;
6760   t->ksp_its           = 0;
6761   t->snes_its          = 0;
6762   t->nwork             = 0;
6763   t->rhsjacobian.time  = -1e20;
6764   t->rhsjacobian.scale = 1.;
6765   t->ijacobian.shift   = 1.;
6766 
6767   ierr = TSGetSNES(tsin,&snes_start);                   CHKERRQ(ierr);
6768   ierr = TSSetSNES(t,snes_start);                       CHKERRQ(ierr);
6769 
6770   ierr = TSGetDM(tsin,&dm);                             CHKERRQ(ierr);
6771   ierr = TSSetDM(t,dm);                                 CHKERRQ(ierr);
6772 
6773   t->adapt=tsin->adapt;
6774   PetscObjectReference((PetscObject)t->adapt);
6775 
6776   t->problem_type      = tsin->problem_type;
6777   t->ptime             = tsin->ptime;
6778   t->time_step         = tsin->time_step;
6779   t->time_step_orig    = tsin->time_step_orig;
6780   t->max_time          = tsin->max_time;
6781   t->steps             = tsin->steps;
6782   t->max_steps         = tsin->max_steps;
6783   t->equation_type     = tsin->equation_type;
6784   t->atol              = tsin->atol;
6785   t->rtol              = tsin->rtol;
6786   t->max_snes_failures = tsin->max_snes_failures;
6787   t->max_reject        = tsin->max_reject;
6788   t->errorifstepfailed = tsin->errorifstepfailed;
6789 
6790   ierr = TSGetType(tsin,&type); CHKERRQ(ierr);
6791   ierr = TSSetType(t,type);     CHKERRQ(ierr);
6792 
6793   t->vec_sol           = NULL;
6794 
6795   t->cfltime          = tsin->cfltime;
6796   t->cfltime_local    = tsin->cfltime_local;
6797   t->exact_final_time = tsin->exact_final_time;
6798 
6799   ierr = PetscMemcpy(t->ops,tsin->ops,sizeof(struct _TSOps));CHKERRQ(ierr);
6800 
6801   if (((PetscObject)tsin)->fortran_func_pointers) {
6802     PetscInt i;
6803     ierr = PetscMalloc((10)*sizeof(void(*)(void)),&((PetscObject)t)->fortran_func_pointers);CHKERRQ(ierr);
6804     for (i=0; i<10; i++) {
6805       ((PetscObject)t)->fortran_func_pointers[i] = ((PetscObject)tsin)->fortran_func_pointers[i];
6806     }
6807   }
6808   *tsout = t;
6809   PetscFunctionReturn(0);
6810 }
6811