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