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