xref: /petsc/src/ts/interface/ts.c (revision 4bb054149a39f68d71e57d43ff252ca13d39c54b)
173f4d377SMatthew Knepley /* $Id: ts.c,v 1.43 2001/09/07 20:12:01 bsmith Exp $ */
2e090d566SSatish Balay #include "src/ts/tsimpl.h"        /*I "petscts.h"  I*/
3d763cef2SBarry Smith 
4d5ba7fb7SMatthew Knepley /* Logging support */
58ba1e511SMatthew Knepley int TS_COOKIE;
6d5ba7fb7SMatthew Knepley int TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
7d405a339SMatthew Knepley 
84a2ae208SSatish Balay #undef __FUNCT__
9bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
10bdad233fSMatthew Knepley /*
11bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
12bdad233fSMatthew Knepley 
13bdad233fSMatthew Knepley   Collective on TS
14bdad233fSMatthew Knepley 
15bdad233fSMatthew Knepley   Input Parameter:
16bdad233fSMatthew Knepley . ts - The ts
17bdad233fSMatthew Knepley 
18bdad233fSMatthew Knepley   Level: intermediate
19bdad233fSMatthew Knepley 
20bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
21bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
22bdad233fSMatthew Knepley */
23bdad233fSMatthew Knepley static int TSSetTypeFromOptions(TS ts)
24bdad233fSMatthew Knepley {
25bdad233fSMatthew Knepley   PetscTruth opt;
26bdad233fSMatthew Knepley   char      *defaultType;
27bdad233fSMatthew Knepley   char       typeName[256];
28bdad233fSMatthew Knepley   int        ierr;
29bdad233fSMatthew Knepley 
30bdad233fSMatthew Knepley   PetscFunctionBegin;
31bdad233fSMatthew Knepley   if (ts->type_name != PETSC_NULL) {
32bdad233fSMatthew Knepley     defaultType = ts->type_name;
33bdad233fSMatthew Knepley   } else {
34bdad233fSMatthew Knepley     defaultType = TS_EULER;
35bdad233fSMatthew Knepley   }
36bdad233fSMatthew Knepley 
37bdad233fSMatthew Knepley   if (!TSRegisterAllCalled) {
38bdad233fSMatthew Knepley     ierr = TSRegisterAll(PETSC_NULL);                                                                     CHKERRQ(ierr);
39bdad233fSMatthew Knepley   }
40bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
41bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
42bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);                                                                       CHKERRQ(ierr);
43bdad233fSMatthew Knepley   } else {
44bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);                                                                    CHKERRQ(ierr);
45bdad233fSMatthew Knepley   }
46bdad233fSMatthew Knepley   PetscFunctionReturn(0);
47bdad233fSMatthew Knepley }
48bdad233fSMatthew Knepley 
49bdad233fSMatthew Knepley #undef __FUNCT__
50bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
51bdad233fSMatthew Knepley /*@
52bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
53bdad233fSMatthew Knepley 
54bdad233fSMatthew Knepley    Collective on TS
55bdad233fSMatthew Knepley 
56bdad233fSMatthew Knepley    Input Parameter:
57bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
58bdad233fSMatthew Knepley 
59bdad233fSMatthew Knepley    Options Database Keys:
60bdad233fSMatthew Knepley +  -ts_type <type> - TS_EULER, TS_BEULER, TS_PVODE, TS_PSEUDO, TS_CRANK_NICHOLSON
61bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
62bdad233fSMatthew Knepley .  -ts_max_time time - maximum time to compute to
63bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
64bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
65bdad233fSMatthew Knepley -  -ts_xmonitor - plot information at each timestep
66bdad233fSMatthew Knepley 
67bdad233fSMatthew Knepley    Level: beginner
68bdad233fSMatthew Knepley 
69bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
70bdad233fSMatthew Knepley 
71bdad233fSMatthew Knepley .seealso: TSGetType
72bdad233fSMatthew Knepley @*/
73bdad233fSMatthew Knepley int TSSetFromOptions(TS ts)
74bdad233fSMatthew Knepley {
75bdad233fSMatthew Knepley   PetscReal  dt;
76bdad233fSMatthew Knepley   PetscTruth opt;
77bdad233fSMatthew Knepley   int        ierr;
78bdad233fSMatthew Knepley 
79bdad233fSMatthew Knepley   PetscFunctionBegin;
80bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
81bdad233fSMatthew Knepley   ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");                              CHKERRQ(ierr);
82bdad233fSMatthew Knepley 
83bdad233fSMatthew Knepley   /* Handle generic TS options */
84bdad233fSMatthew Knepley   ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
85bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
86bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
87bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
88bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
89bdad233fSMatthew Knepley     ts->initial_time_step = ts->time_step = dt;
90bdad233fSMatthew Knepley   }
91bdad233fSMatthew Knepley 
92bdad233fSMatthew Knepley   /* Monitor options */
93bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_monitor","Monitor timestep size","TSDefaultMonitor",&opt);               CHKERRQ(ierr);
94bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
95bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSDefaultMonitor,PETSC_NULL,PETSC_NULL);                                     CHKERRQ(ierr);
96bdad233fSMatthew Knepley     }
97bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_xmonitor","Monitor timestep size graphically","TSLGMonitor",&opt);       CHKERRQ(ierr);
98bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
99bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSLGMonitor,PETSC_NULL,PETSC_NULL);                                          CHKERRQ(ierr);
100bdad233fSMatthew Knepley     }
101bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_vecmonitor","Monitor solution graphically","TSVecViewMonitor",&opt);     CHKERRQ(ierr);
102bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
103bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSVecViewMonitor,PETSC_NULL,PETSC_NULL);                                     CHKERRQ(ierr);
104bdad233fSMatthew Knepley     }
105bdad233fSMatthew Knepley 
106bdad233fSMatthew Knepley   /* Handle TS type options */
107bdad233fSMatthew Knepley   ierr = TSSetTypeFromOptions(ts);                                                                        CHKERRQ(ierr);
108bdad233fSMatthew Knepley 
109bdad233fSMatthew Knepley   /* Handle specific TS options */
110bdad233fSMatthew Knepley   if (ts->ops->setfromoptions != PETSC_NULL) {
111bdad233fSMatthew Knepley     ierr = (*ts->ops->setfromoptions)(ts);                                                                CHKERRQ(ierr);
112bdad233fSMatthew Knepley   }
113bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();                                                                               CHKERRQ(ierr);
114bdad233fSMatthew Knepley 
115bdad233fSMatthew Knepley   /* Handle subobject options */
116bdad233fSMatthew Knepley   switch(ts->problem_type) {
117156fc9a6SMatthew Knepley     /* Should check for implicit/explicit */
118bdad233fSMatthew Knepley   case TS_LINEAR:
119156fc9a6SMatthew Knepley     if (ts->sles != PETSC_NULL) {
120bdad233fSMatthew Knepley       ierr = SLESSetFromOptions(ts->sles);                                                                CHKERRQ(ierr);
121156fc9a6SMatthew Knepley     }
122bdad233fSMatthew Knepley     break;
123bdad233fSMatthew Knepley   case TS_NONLINEAR:
124156fc9a6SMatthew Knepley     if (ts->snes != PETSC_NULL) {
125bdad233fSMatthew Knepley       ierr = SNESSetFromOptions(ts->snes);                                                                CHKERRQ(ierr);
126156fc9a6SMatthew Knepley     }
127bdad233fSMatthew Knepley     break;
128bdad233fSMatthew Knepley   default:
129bdad233fSMatthew Knepley     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", ts->problem_type);
130bdad233fSMatthew Knepley   }
131bdad233fSMatthew Knepley 
132bdad233fSMatthew Knepley   PetscFunctionReturn(0);
133bdad233fSMatthew Knepley }
134bdad233fSMatthew Knepley 
135bdad233fSMatthew Knepley #undef  __FUNCT__
136bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
137bdad233fSMatthew Knepley /*@
138bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
139bdad233fSMatthew Knepley 
140bdad233fSMatthew Knepley   Collective on TS
141bdad233fSMatthew Knepley 
142bdad233fSMatthew Knepley   Input Parameter:
143bdad233fSMatthew Knepley . ts - The ts
144bdad233fSMatthew Knepley 
145bdad233fSMatthew Knepley   Level: intermediate
146bdad233fSMatthew Knepley 
147bdad233fSMatthew Knepley .keywords: TS, view, options, database
148bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
149bdad233fSMatthew Knepley @*/
150bdad233fSMatthew Knepley int TSViewFromOptions(TS ts, char *title)
151bdad233fSMatthew Knepley {
152bdad233fSMatthew Knepley   PetscViewer viewer;
153bdad233fSMatthew Knepley   PetscDraw   draw;
154bdad233fSMatthew Knepley   PetscTruth  opt;
155bdad233fSMatthew Knepley   char       *titleStr;
156bdad233fSMatthew Knepley   char        typeName[1024];
157bdad233fSMatthew Knepley   char        fileName[1024];
158bdad233fSMatthew Knepley   int         len;
159bdad233fSMatthew Knepley   int         ierr;
160bdad233fSMatthew Knepley 
161bdad233fSMatthew Knepley   PetscFunctionBegin;
162bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);                                               CHKERRQ(ierr);
163bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
164bdad233fSMatthew Knepley     ierr = PetscOptionsGetString(ts->prefix, "-ts_view", typeName, 1024, &opt);                           CHKERRQ(ierr);
165bdad233fSMatthew Knepley     ierr = PetscStrlen(typeName, &len);                                                                   CHKERRQ(ierr);
166bdad233fSMatthew Knepley     if (len > 0) {
167bdad233fSMatthew Knepley       ierr = PetscViewerCreate(ts->comm, &viewer);                                                        CHKERRQ(ierr);
168bdad233fSMatthew Knepley       ierr = PetscViewerSetType(viewer, typeName);                                                        CHKERRQ(ierr);
169bdad233fSMatthew Knepley       ierr = PetscOptionsGetString(ts->prefix, "-ts_view_file", fileName, 1024, &opt);                    CHKERRQ(ierr);
170bdad233fSMatthew Knepley       if (opt == PETSC_TRUE) {
171bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, fileName);                                                  CHKERRQ(ierr);
172bdad233fSMatthew Knepley       } else {
173bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, ts->name);                                                  CHKERRQ(ierr);
174bdad233fSMatthew Knepley       }
175bdad233fSMatthew Knepley       ierr = TSView(ts, viewer);                                                                          CHKERRQ(ierr);
176bdad233fSMatthew Knepley       ierr = PetscViewerFlush(viewer);                                                                    CHKERRQ(ierr);
177bdad233fSMatthew Knepley       ierr = PetscViewerDestroy(viewer);                                                                  CHKERRQ(ierr);
178bdad233fSMatthew Knepley     } else {
179bdad233fSMatthew Knepley       ierr = TSView(ts, PETSC_NULL);                                                                      CHKERRQ(ierr);
180bdad233fSMatthew Knepley     }
181bdad233fSMatthew Knepley   }
182bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);                                          CHKERRQ(ierr);
183bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
184bdad233fSMatthew Knepley     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);                                  CHKERRQ(ierr);
185bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);                                                      CHKERRQ(ierr);
186bdad233fSMatthew Knepley     if (title != PETSC_NULL) {
187bdad233fSMatthew Knepley       titleStr = title;
188bdad233fSMatthew Knepley     } else {
189bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject) ts);                                                           CHKERRQ(ierr) ;
190bdad233fSMatthew Knepley       titleStr = ts->name;
191bdad233fSMatthew Knepley     }
192bdad233fSMatthew Knepley     ierr = PetscDrawSetTitle(draw, titleStr);                                                             CHKERRQ(ierr);
193bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);                                                                            CHKERRQ(ierr);
194bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);                                                                      CHKERRQ(ierr);
195bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);                                                                          CHKERRQ(ierr);
196bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);                                                                    CHKERRQ(ierr);
197bdad233fSMatthew Knepley   }
198bdad233fSMatthew Knepley   PetscFunctionReturn(0);
199bdad233fSMatthew Knepley }
200bdad233fSMatthew Knepley 
201bdad233fSMatthew Knepley #undef __FUNCT__
2024a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
203a7a1495cSBarry Smith /*@
2048c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
205a7a1495cSBarry Smith       set with TSSetRHSJacobian().
206a7a1495cSBarry Smith 
207a7a1495cSBarry Smith    Collective on TS and Vec
208a7a1495cSBarry Smith 
209a7a1495cSBarry Smith    Input Parameters:
210a7a1495cSBarry Smith +  ts - the SNES context
211a7a1495cSBarry Smith .  t - current timestep
212a7a1495cSBarry Smith -  x - input vector
213a7a1495cSBarry Smith 
214a7a1495cSBarry Smith    Output Parameters:
215a7a1495cSBarry Smith +  A - Jacobian matrix
216a7a1495cSBarry Smith .  B - optional preconditioning matrix
217a7a1495cSBarry Smith -  flag - flag indicating matrix structure
218a7a1495cSBarry Smith 
219a7a1495cSBarry Smith    Notes:
220a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
221a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
222a7a1495cSBarry Smith 
223a7a1495cSBarry Smith    See SLESSetOperators() for important information about setting the
224a7a1495cSBarry Smith    flag parameter.
225a7a1495cSBarry Smith 
226a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
227a7a1495cSBarry Smith 
228a7a1495cSBarry Smith    Level: developer
229a7a1495cSBarry Smith 
230a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
231a7a1495cSBarry Smith 
232a7a1495cSBarry Smith .seealso:  TSSetRHSJacobian(), SLESSetOperators()
233a7a1495cSBarry Smith @*/
23487828ca2SBarry Smith int TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
235a7a1495cSBarry Smith {
236a7a1495cSBarry Smith   int ierr;
237a7a1495cSBarry Smith 
238a7a1495cSBarry Smith   PetscFunctionBegin;
239a7a1495cSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
240a7a1495cSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE);
241a7a1495cSBarry Smith   PetscCheckSameComm(ts,X);
242a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
24329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
244a7a1495cSBarry Smith   }
245000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
246d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
247a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
248a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
249000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
250a7a1495cSBarry Smith     PetscStackPop;
251d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
252a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
253a7a1495cSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE);
254a7a1495cSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE);
255ef66eb69SBarry Smith   } else {
256ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258ef66eb69SBarry Smith     if (*A != *B) {
259ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
260ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261ef66eb69SBarry Smith     }
262ef66eb69SBarry Smith   }
263a7a1495cSBarry Smith   PetscFunctionReturn(0);
264a7a1495cSBarry Smith }
265a7a1495cSBarry Smith 
2664a2ae208SSatish Balay #undef __FUNCT__
2674a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
268d763cef2SBarry Smith /*
269d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
270d763cef2SBarry Smith 
271d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
272d763cef2SBarry Smith    this routine applies the matrix.
273d763cef2SBarry Smith */
27487828ca2SBarry Smith int TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
275d763cef2SBarry Smith {
276d763cef2SBarry Smith   int ierr;
277d763cef2SBarry Smith 
278d763cef2SBarry Smith   PetscFunctionBegin;
279d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
2807c922b88SBarry Smith   PetscValidHeader(x);
2817c922b88SBarry Smith   PetscValidHeader(y);
282d763cef2SBarry Smith 
283d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
284000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
285d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
286000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
287d763cef2SBarry Smith     PetscStackPop;
2887c922b88SBarry Smith   } else {
289000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
290d763cef2SBarry Smith       MatStructure flg;
291d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
292000e7ae3SMatthew Knepley       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
293d763cef2SBarry Smith       PetscStackPop;
294d763cef2SBarry Smith     }
295d763cef2SBarry Smith     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
2967c922b88SBarry Smith   }
297d763cef2SBarry Smith 
298d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
299d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
300d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
301d763cef2SBarry Smith 
302d763cef2SBarry Smith   PetscFunctionReturn(0);
303d763cef2SBarry Smith }
304d763cef2SBarry Smith 
3054a2ae208SSatish Balay #undef __FUNCT__
3064a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
307d763cef2SBarry Smith /*@C
308d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
309d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
310d763cef2SBarry Smith 
311d763cef2SBarry Smith     Collective on TS
312d763cef2SBarry Smith 
313d763cef2SBarry Smith     Input Parameters:
314d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
315d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
316d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
317d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
318d763cef2SBarry Smith 
319d763cef2SBarry Smith     Calling sequence of func:
32087828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
321d763cef2SBarry Smith 
322d763cef2SBarry Smith +   t - current timestep
323d763cef2SBarry Smith .   u - input vector
324d763cef2SBarry Smith .   F - function vector
325d763cef2SBarry Smith -   ctx - [optional] user-defined function context
326d763cef2SBarry Smith 
327d763cef2SBarry Smith     Important:
328d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
329d763cef2SBarry Smith 
330d763cef2SBarry Smith     Level: beginner
331d763cef2SBarry Smith 
332d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
333d763cef2SBarry Smith 
334d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
335d763cef2SBarry Smith @*/
33687828ca2SBarry Smith int TSSetRHSFunction(TS ts,int (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
337d763cef2SBarry Smith {
338d763cef2SBarry Smith   PetscFunctionBegin;
339d763cef2SBarry Smith 
340d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
341d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
34229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
343d763cef2SBarry Smith   }
344000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
345d763cef2SBarry Smith   ts->funP             = ctx;
346d763cef2SBarry Smith   PetscFunctionReturn(0);
347d763cef2SBarry Smith }
348d763cef2SBarry Smith 
3494a2ae208SSatish Balay #undef __FUNCT__
3504a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
351d763cef2SBarry Smith /*@C
352d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
353d763cef2SBarry Smith    Also sets the location to store A.
354d763cef2SBarry Smith 
355d763cef2SBarry Smith    Collective on TS
356d763cef2SBarry Smith 
357d763cef2SBarry Smith    Input Parameters:
358d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
359d763cef2SBarry Smith .  A   - matrix
360d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
361d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
362d763cef2SBarry Smith          if A is not a function of t.
363d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
364d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
365d763cef2SBarry Smith 
366d763cef2SBarry Smith    Calling sequence of func:
36787828ca2SBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,int *flag,void *ctx);
368d763cef2SBarry Smith 
369d763cef2SBarry Smith +  t - current timestep
370d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
371d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
372d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
373d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
374d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
375d763cef2SBarry Smith 
376d763cef2SBarry Smith    Notes:
377d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
378d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
379d763cef2SBarry Smith 
380d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
381d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
382d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
383d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
384d763cef2SBarry Smith    throughout the global iterations.
385d763cef2SBarry Smith 
386d763cef2SBarry Smith    Important:
387d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
388d763cef2SBarry Smith 
389d763cef2SBarry Smith    Level: beginner
390d763cef2SBarry Smith 
391d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
392d763cef2SBarry Smith 
393d763cef2SBarry Smith .seealso: TSSetRHSFunction()
394d763cef2SBarry Smith @*/
39587828ca2SBarry Smith int TSSetRHSMatrix(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
396d763cef2SBarry Smith {
397d763cef2SBarry Smith   PetscFunctionBegin;
398d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
399184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
400184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
401184914b5SBarry Smith   PetscCheckSameComm(ts,A);
402184914b5SBarry Smith   PetscCheckSameComm(ts,B);
403d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
40429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
405d763cef2SBarry Smith   }
406d763cef2SBarry Smith 
407000e7ae3SMatthew Knepley   ts->ops->rhsmatrix = f;
408d763cef2SBarry Smith   ts->jacP           = ctx;
409d763cef2SBarry Smith   ts->A              = A;
410d763cef2SBarry Smith   ts->B              = B;
411d763cef2SBarry Smith 
412d763cef2SBarry Smith   PetscFunctionReturn(0);
413d763cef2SBarry Smith }
414d763cef2SBarry Smith 
4154a2ae208SSatish Balay #undef __FUNCT__
4164a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
417d763cef2SBarry Smith /*@C
418d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
419d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
420d763cef2SBarry Smith 
421d763cef2SBarry Smith    Collective on TS
422d763cef2SBarry Smith 
423d763cef2SBarry Smith    Input Parameters:
424d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
425d763cef2SBarry Smith .  A   - Jacobian matrix
426d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
427d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
428d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
429d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
430d763cef2SBarry Smith 
431d763cef2SBarry Smith    Calling sequence of func:
43287828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
433d763cef2SBarry Smith 
434d763cef2SBarry Smith +  t - current timestep
435d763cef2SBarry Smith .  u - input vector
436d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
437d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
438d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
439d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
440d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
441d763cef2SBarry Smith 
442d763cef2SBarry Smith    Notes:
443d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
444d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
445d763cef2SBarry Smith 
446d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
447d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
448d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
449d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
450d763cef2SBarry Smith    throughout the global iterations.
451d763cef2SBarry Smith 
452d763cef2SBarry Smith    Level: beginner
453d763cef2SBarry Smith 
454d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
455d763cef2SBarry Smith 
456d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
457d763cef2SBarry Smith           SNESDefaultComputeJacobianColor()
458d763cef2SBarry Smith 
459d763cef2SBarry Smith @*/
46087828ca2SBarry Smith int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
461d763cef2SBarry Smith {
462d763cef2SBarry Smith   PetscFunctionBegin;
463d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
464184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
465184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
466184914b5SBarry Smith   PetscCheckSameComm(ts,A);
467184914b5SBarry Smith   PetscCheckSameComm(ts,B);
468d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
46929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
470d763cef2SBarry Smith   }
471d763cef2SBarry Smith 
472000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
473d763cef2SBarry Smith   ts->jacP             = ctx;
474d763cef2SBarry Smith   ts->A                = A;
475d763cef2SBarry Smith   ts->B                = B;
476d763cef2SBarry Smith   PetscFunctionReturn(0);
477d763cef2SBarry Smith }
478d763cef2SBarry Smith 
4794a2ae208SSatish Balay #undef __FUNCT__
4804a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
481d763cef2SBarry Smith /*
482d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
483d763cef2SBarry Smith 
484d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
485d763cef2SBarry Smith    this routine applies the matrix.
486d763cef2SBarry Smith */
48787828ca2SBarry Smith int TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
488d763cef2SBarry Smith {
489d763cef2SBarry Smith   int ierr;
490d763cef2SBarry Smith 
491d763cef2SBarry Smith   PetscFunctionBegin;
492d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
493d763cef2SBarry Smith   PetscValidHeader(x);
494184914b5SBarry Smith   PetscCheckSameComm(ts,x);
495d763cef2SBarry Smith 
496000e7ae3SMatthew Knepley   if (ts->ops->rhsbc) {
497d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
498000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
499d763cef2SBarry Smith     PetscStackPop;
500d763cef2SBarry Smith     PetscFunctionReturn(0);
501d763cef2SBarry Smith   }
502d763cef2SBarry Smith 
503d763cef2SBarry Smith   PetscFunctionReturn(0);
504d763cef2SBarry Smith }
505d763cef2SBarry Smith 
5064a2ae208SSatish Balay #undef __FUNCT__
5074a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
508d763cef2SBarry Smith /*@C
509d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
510d763cef2SBarry Smith     boundary conditions for the function F.
511d763cef2SBarry Smith 
512d763cef2SBarry Smith     Collective on TS
513d763cef2SBarry Smith 
514d763cef2SBarry Smith     Input Parameters:
515d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
516d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
517d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
518d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
519d763cef2SBarry Smith 
520d763cef2SBarry Smith     Calling sequence of func:
52187828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
522d763cef2SBarry Smith 
523d763cef2SBarry Smith +   t - current timestep
524d763cef2SBarry Smith .   F - function vector
525d763cef2SBarry Smith -   ctx - [optional] user-defined function context
526d763cef2SBarry Smith 
527d763cef2SBarry Smith     Level: intermediate
528d763cef2SBarry Smith 
529d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
530d763cef2SBarry Smith @*/
53187828ca2SBarry Smith int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,PetscReal,Vec,void*),void *ctx)
532d763cef2SBarry Smith {
533d763cef2SBarry Smith   PetscFunctionBegin;
534d763cef2SBarry Smith 
535d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
536d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
53729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
538d763cef2SBarry Smith   }
539000e7ae3SMatthew Knepley   ts->ops->rhsbc = f;
540d763cef2SBarry Smith   ts->bcP        = ctx;
541d763cef2SBarry Smith   PetscFunctionReturn(0);
542d763cef2SBarry Smith }
543d763cef2SBarry Smith 
5444a2ae208SSatish Balay #undef __FUNCT__
5454a2ae208SSatish Balay #define __FUNCT__ "TSView"
5467e2c5f70SBarry Smith /*@C
547d763cef2SBarry Smith     TSView - Prints the TS data structure.
548d763cef2SBarry Smith 
5494c49b128SBarry Smith     Collective on TS
550d763cef2SBarry Smith 
551d763cef2SBarry Smith     Input Parameters:
552d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
553d763cef2SBarry Smith -   viewer - visualization context
554d763cef2SBarry Smith 
555d763cef2SBarry Smith     Options Database Key:
556d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
557d763cef2SBarry Smith 
558d763cef2SBarry Smith     Notes:
559d763cef2SBarry Smith     The available visualization contexts include
560b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
561b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
562d763cef2SBarry Smith          output where only the first processor opens
563d763cef2SBarry Smith          the file.  All other processors send their
564d763cef2SBarry Smith          data to the first processor to print.
565d763cef2SBarry Smith 
566d763cef2SBarry Smith     The user can open an alternative visualization context with
567b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
568d763cef2SBarry Smith 
569d763cef2SBarry Smith     Level: beginner
570d763cef2SBarry Smith 
571d763cef2SBarry Smith .keywords: TS, timestep, view
572d763cef2SBarry Smith 
573b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
574d763cef2SBarry Smith @*/
575b0a32e0cSBarry Smith int TSView(TS ts,PetscViewer viewer)
576d763cef2SBarry Smith {
577d763cef2SBarry Smith   int        ierr;
578454a90a3SBarry Smith   char       *type;
5796831982aSBarry Smith   PetscTruth isascii,isstring;
580d763cef2SBarry Smith 
581d763cef2SBarry Smith   PetscFunctionBegin;
582d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
583b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
584b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
5856831982aSBarry Smith   PetscCheckSameComm(ts,viewer);
586fd16b177SBarry Smith 
587b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
588b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
5890f5bd95cSBarry Smith   if (isascii) {
590b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
591454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
592454a90a3SBarry Smith     if (type) {
593b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
594184914b5SBarry Smith     } else {
595b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
596184914b5SBarry Smith     }
597000e7ae3SMatthew Knepley     if (ts->ops->view) {
598b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
599000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
600b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
601d763cef2SBarry Smith     }
602b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
603b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
604d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
605b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
606d763cef2SBarry Smith     }
607b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
6080f5bd95cSBarry Smith   } else if (isstring) {
609454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
610b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
611d763cef2SBarry Smith   }
612b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
613d763cef2SBarry Smith   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
614d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
615b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
616d763cef2SBarry Smith   PetscFunctionReturn(0);
617d763cef2SBarry Smith }
618d763cef2SBarry Smith 
619d763cef2SBarry Smith 
6204a2ae208SSatish Balay #undef __FUNCT__
6214a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
622d763cef2SBarry Smith /*@C
623d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
624d763cef2SBarry Smith    the timesteppers.
625d763cef2SBarry Smith 
626d763cef2SBarry Smith    Collective on TS
627d763cef2SBarry Smith 
628d763cef2SBarry Smith    Input Parameters:
629d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
630d763cef2SBarry Smith -  usrP - optional user context
631d763cef2SBarry Smith 
632d763cef2SBarry Smith    Level: intermediate
633d763cef2SBarry Smith 
634d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
635d763cef2SBarry Smith 
636d763cef2SBarry Smith .seealso: TSGetApplicationContext()
637d763cef2SBarry Smith @*/
638d763cef2SBarry Smith int TSSetApplicationContext(TS ts,void *usrP)
639d763cef2SBarry Smith {
640d763cef2SBarry Smith   PetscFunctionBegin;
641d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
642d763cef2SBarry Smith   ts->user = usrP;
643d763cef2SBarry Smith   PetscFunctionReturn(0);
644d763cef2SBarry Smith }
645d763cef2SBarry Smith 
6464a2ae208SSatish Balay #undef __FUNCT__
6474a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
648d763cef2SBarry Smith /*@C
649d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
650d763cef2SBarry Smith     timestepper.
651d763cef2SBarry Smith 
652d763cef2SBarry Smith     Not Collective
653d763cef2SBarry Smith 
654d763cef2SBarry Smith     Input Parameter:
655d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
656d763cef2SBarry Smith 
657d763cef2SBarry Smith     Output Parameter:
658d763cef2SBarry Smith .   usrP - user context
659d763cef2SBarry Smith 
660d763cef2SBarry Smith     Level: intermediate
661d763cef2SBarry Smith 
662d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
663d763cef2SBarry Smith 
664d763cef2SBarry Smith .seealso: TSSetApplicationContext()
665d763cef2SBarry Smith @*/
666d763cef2SBarry Smith int TSGetApplicationContext(TS ts,void **usrP)
667d763cef2SBarry Smith {
668d763cef2SBarry Smith   PetscFunctionBegin;
669d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
670d763cef2SBarry Smith   *usrP = ts->user;
671d763cef2SBarry Smith   PetscFunctionReturn(0);
672d763cef2SBarry Smith }
673d763cef2SBarry Smith 
6744a2ae208SSatish Balay #undef __FUNCT__
6754a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
676d763cef2SBarry Smith /*@
677d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
678d763cef2SBarry Smith 
679d763cef2SBarry Smith    Not Collective
680d763cef2SBarry Smith 
681d763cef2SBarry Smith    Input Parameter:
682d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
683d763cef2SBarry Smith 
684d763cef2SBarry Smith    Output Parameter:
685d763cef2SBarry Smith .  iter - number steps so far
686d763cef2SBarry Smith 
687d763cef2SBarry Smith    Level: intermediate
688d763cef2SBarry Smith 
689d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
690d763cef2SBarry Smith @*/
691d763cef2SBarry Smith int TSGetTimeStepNumber(TS ts,int* iter)
692d763cef2SBarry Smith {
693d763cef2SBarry Smith   PetscFunctionBegin;
694d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
695d763cef2SBarry Smith   *iter = ts->steps;
696d763cef2SBarry Smith   PetscFunctionReturn(0);
697d763cef2SBarry Smith }
698d763cef2SBarry Smith 
6994a2ae208SSatish Balay #undef __FUNCT__
7004a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
701d763cef2SBarry Smith /*@
702d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
703d763cef2SBarry Smith    as well as the initial time.
704d763cef2SBarry Smith 
705d763cef2SBarry Smith    Collective on TS
706d763cef2SBarry Smith 
707d763cef2SBarry Smith    Input Parameters:
708d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
709d763cef2SBarry Smith .  initial_time - the initial time
710d763cef2SBarry Smith -  time_step - the size of the timestep
711d763cef2SBarry Smith 
712d763cef2SBarry Smith    Level: intermediate
713d763cef2SBarry Smith 
714d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
715d763cef2SBarry Smith 
716d763cef2SBarry Smith .keywords: TS, set, initial, timestep
717d763cef2SBarry Smith @*/
71887828ca2SBarry Smith int TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
719d763cef2SBarry Smith {
720d763cef2SBarry Smith   PetscFunctionBegin;
721d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
722d763cef2SBarry Smith   ts->time_step         = time_step;
723d763cef2SBarry Smith   ts->initial_time_step = time_step;
724d763cef2SBarry Smith   ts->ptime             = initial_time;
725d763cef2SBarry Smith   PetscFunctionReturn(0);
726d763cef2SBarry Smith }
727d763cef2SBarry Smith 
7284a2ae208SSatish Balay #undef __FUNCT__
7294a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
730d763cef2SBarry Smith /*@
731d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
732d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
733d763cef2SBarry Smith 
734d763cef2SBarry Smith    Collective on TS
735d763cef2SBarry Smith 
736d763cef2SBarry Smith    Input Parameters:
737d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
738d763cef2SBarry Smith -  time_step - the size of the timestep
739d763cef2SBarry Smith 
740d763cef2SBarry Smith    Level: intermediate
741d763cef2SBarry Smith 
742d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
743d763cef2SBarry Smith 
744d763cef2SBarry Smith .keywords: TS, set, timestep
745d763cef2SBarry Smith @*/
74687828ca2SBarry Smith int TSSetTimeStep(TS ts,PetscReal time_step)
747d763cef2SBarry Smith {
748d763cef2SBarry Smith   PetscFunctionBegin;
749d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
750d763cef2SBarry Smith   ts->time_step = time_step;
751d763cef2SBarry Smith   PetscFunctionReturn(0);
752d763cef2SBarry Smith }
753d763cef2SBarry Smith 
7544a2ae208SSatish Balay #undef __FUNCT__
7554a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
756d763cef2SBarry Smith /*@
757d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
758d763cef2SBarry Smith 
759d763cef2SBarry Smith    Not Collective
760d763cef2SBarry Smith 
761d763cef2SBarry Smith    Input Parameter:
762d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
763d763cef2SBarry Smith 
764d763cef2SBarry Smith    Output Parameter:
765d763cef2SBarry Smith .  dt - the current timestep size
766d763cef2SBarry Smith 
767d763cef2SBarry Smith    Level: intermediate
768d763cef2SBarry Smith 
769d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
770d763cef2SBarry Smith 
771d763cef2SBarry Smith .keywords: TS, get, timestep
772d763cef2SBarry Smith @*/
77387828ca2SBarry Smith int TSGetTimeStep(TS ts,PetscReal* dt)
774d763cef2SBarry Smith {
775d763cef2SBarry Smith   PetscFunctionBegin;
776d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
777d763cef2SBarry Smith   *dt = ts->time_step;
778d763cef2SBarry Smith   PetscFunctionReturn(0);
779d763cef2SBarry Smith }
780d763cef2SBarry Smith 
7814a2ae208SSatish Balay #undef __FUNCT__
7824a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
783d763cef2SBarry Smith /*@C
784d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
785d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
786d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
787d763cef2SBarry Smith    the solution at the next timestep has been calculated.
788d763cef2SBarry Smith 
789d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
790d763cef2SBarry Smith 
791d763cef2SBarry Smith    Input Parameter:
792d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
793d763cef2SBarry Smith 
794d763cef2SBarry Smith    Output Parameter:
795d763cef2SBarry Smith .  v - the vector containing the solution
796d763cef2SBarry Smith 
797d763cef2SBarry Smith    Level: intermediate
798d763cef2SBarry Smith 
799d763cef2SBarry Smith .seealso: TSGetTimeStep()
800d763cef2SBarry Smith 
801d763cef2SBarry Smith .keywords: TS, timestep, get, solution
802d763cef2SBarry Smith @*/
803d763cef2SBarry Smith int TSGetSolution(TS ts,Vec *v)
804d763cef2SBarry Smith {
805d763cef2SBarry Smith   PetscFunctionBegin;
806d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
807d763cef2SBarry Smith   *v = ts->vec_sol_always;
808d763cef2SBarry Smith   PetscFunctionReturn(0);
809d763cef2SBarry Smith }
810d763cef2SBarry Smith 
811bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
8124a2ae208SSatish Balay #undef __FUNCT__
813bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
814d763cef2SBarry Smith /*@C
815bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
816d763cef2SBarry Smith 
817bdad233fSMatthew Knepley   Not collective
818d763cef2SBarry Smith 
819bdad233fSMatthew Knepley   Input Parameters:
820bdad233fSMatthew Knepley + ts   - The TS
821bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
822d763cef2SBarry Smith .vb
823d763cef2SBarry Smith          U_t = A U
824d763cef2SBarry Smith          U_t = A(t) U
825d763cef2SBarry Smith          U_t = F(t,U)
826d763cef2SBarry Smith .ve
827d763cef2SBarry Smith 
828d763cef2SBarry Smith    Level: beginner
829d763cef2SBarry Smith 
830bdad233fSMatthew Knepley .keywords: TS, problem type
831bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
832d763cef2SBarry Smith @*/
833bdad233fSMatthew Knepley int TSSetProblemType(TS ts, TSProblemType type) {
834d763cef2SBarry Smith   PetscFunctionBegin;
835bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
836bdad233fSMatthew Knepley   ts->problem_type = type;
837d763cef2SBarry Smith   PetscFunctionReturn(0);
838d763cef2SBarry Smith }
839d763cef2SBarry Smith 
840bdad233fSMatthew Knepley #undef __FUNCT__
841bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
842bdad233fSMatthew Knepley /*@C
843bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
844bdad233fSMatthew Knepley 
845bdad233fSMatthew Knepley   Not collective
846bdad233fSMatthew Knepley 
847bdad233fSMatthew Knepley   Input Parameter:
848bdad233fSMatthew Knepley . ts   - The TS
849bdad233fSMatthew Knepley 
850bdad233fSMatthew Knepley   Output Parameter:
851bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
852bdad233fSMatthew Knepley .vb
853bdad233fSMatthew Knepley          U_t = A U
854bdad233fSMatthew Knepley          U_t = A(t) U
855bdad233fSMatthew Knepley          U_t = F(t,U)
856bdad233fSMatthew Knepley .ve
857bdad233fSMatthew Knepley 
858bdad233fSMatthew Knepley    Level: beginner
859bdad233fSMatthew Knepley 
860bdad233fSMatthew Knepley .keywords: TS, problem type
861bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
862bdad233fSMatthew Knepley @*/
863bdad233fSMatthew Knepley int TSGetProblemType(TS ts, TSProblemType *type) {
864bdad233fSMatthew Knepley   PetscFunctionBegin;
865bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
866bdad233fSMatthew Knepley   PetscValidPointer(type);
867bdad233fSMatthew Knepley   *type = ts->problem_type;
868bdad233fSMatthew Knepley   PetscFunctionReturn(0);
869bdad233fSMatthew Knepley }
870d763cef2SBarry Smith 
8714a2ae208SSatish Balay #undef __FUNCT__
8724a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
873d763cef2SBarry Smith /*@
874d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
875d763cef2SBarry Smith    of a timestepper.
876d763cef2SBarry Smith 
877d763cef2SBarry Smith    Collective on TS
878d763cef2SBarry Smith 
879d763cef2SBarry Smith    Input Parameter:
880d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
881d763cef2SBarry Smith 
882d763cef2SBarry Smith    Notes:
883d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
884d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
885d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
886d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
887d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
888d763cef2SBarry Smith 
889d763cef2SBarry Smith    Level: advanced
890d763cef2SBarry Smith 
891d763cef2SBarry Smith .keywords: TS, timestep, setup
892d763cef2SBarry Smith 
893d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
894d763cef2SBarry Smith @*/
895d763cef2SBarry Smith int TSSetUp(TS ts)
896d763cef2SBarry Smith {
897d763cef2SBarry Smith   int ierr;
898d763cef2SBarry Smith 
899d763cef2SBarry Smith   PetscFunctionBegin;
900d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
90129bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
902d763cef2SBarry Smith   if (!ts->type_name) {
903d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
904d763cef2SBarry Smith   }
905000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
906d763cef2SBarry Smith   ts->setupcalled = 1;
907d763cef2SBarry Smith   PetscFunctionReturn(0);
908d763cef2SBarry Smith }
909d763cef2SBarry Smith 
9104a2ae208SSatish Balay #undef __FUNCT__
9114a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
912d763cef2SBarry Smith /*@C
913d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
914d763cef2SBarry Smith    with TSCreate().
915d763cef2SBarry Smith 
916d763cef2SBarry Smith    Collective on TS
917d763cef2SBarry Smith 
918d763cef2SBarry Smith    Input Parameter:
919d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
920d763cef2SBarry Smith 
921d763cef2SBarry Smith    Level: beginner
922d763cef2SBarry Smith 
923d763cef2SBarry Smith .keywords: TS, timestepper, destroy
924d763cef2SBarry Smith 
925d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
926d763cef2SBarry Smith @*/
927d763cef2SBarry Smith int TSDestroy(TS ts)
928d763cef2SBarry Smith {
929329f5518SBarry Smith   int ierr,i;
930d763cef2SBarry Smith 
931d763cef2SBarry Smith   PetscFunctionBegin;
932d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
933d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
934d763cef2SBarry Smith 
935be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
9360f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
937be0abb6dSBarry Smith 
938d763cef2SBarry Smith   if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);}
939d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
940000e7ae3SMatthew Knepley   ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);
941329f5518SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
942329f5518SBarry Smith     if (ts->mdestroy[i]) {
943329f5518SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
944329f5518SBarry Smith     }
945329f5518SBarry Smith   }
946b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)ts);
947d763cef2SBarry Smith   PetscHeaderDestroy((PetscObject)ts);
948d763cef2SBarry Smith   PetscFunctionReturn(0);
949d763cef2SBarry Smith }
950d763cef2SBarry Smith 
9514a2ae208SSatish Balay #undef __FUNCT__
9524a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
953d763cef2SBarry Smith /*@C
954d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
955d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
956d763cef2SBarry Smith 
957d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
958d763cef2SBarry Smith 
959d763cef2SBarry Smith    Input Parameter:
960d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
961d763cef2SBarry Smith 
962d763cef2SBarry Smith    Output Parameter:
963d763cef2SBarry Smith .  snes - the nonlinear solver context
964d763cef2SBarry Smith 
965d763cef2SBarry Smith    Notes:
966d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
967d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
968d763cef2SBarry Smith    SLES, KSP, and PC contexts as well.
969d763cef2SBarry Smith 
970d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
971d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
972d763cef2SBarry Smith 
973d763cef2SBarry Smith    Level: beginner
974d763cef2SBarry Smith 
975d763cef2SBarry Smith .keywords: timestep, get, SNES
976d763cef2SBarry Smith @*/
977d763cef2SBarry Smith int TSGetSNES(TS ts,SNES *snes)
978d763cef2SBarry Smith {
979d763cef2SBarry Smith   PetscFunctionBegin;
980d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
98129bbc08cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetSLES()");
982d763cef2SBarry Smith   *snes = ts->snes;
983d763cef2SBarry Smith   PetscFunctionReturn(0);
984d763cef2SBarry Smith }
985d763cef2SBarry Smith 
9864a2ae208SSatish Balay #undef __FUNCT__
9874a2ae208SSatish Balay #define __FUNCT__ "TSGetSLES"
988d763cef2SBarry Smith /*@C
989d763cef2SBarry Smith    TSGetSLES - Returns the SLES (linear solver) associated with
990d763cef2SBarry Smith    a TS (timestepper) context.
991d763cef2SBarry Smith 
992d763cef2SBarry Smith    Not Collective, but SLES is parallel if TS is parallel
993d763cef2SBarry Smith 
994d763cef2SBarry Smith    Input Parameter:
995d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
996d763cef2SBarry Smith 
997d763cef2SBarry Smith    Output Parameter:
998d763cef2SBarry Smith .  sles - the nonlinear solver context
999d763cef2SBarry Smith 
1000d763cef2SBarry Smith    Notes:
1001d763cef2SBarry Smith    The user can then directly manipulate the SLES context to set various
1002d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1003d763cef2SBarry Smith    KSP and PC contexts as well.
1004d763cef2SBarry Smith 
1005d763cef2SBarry Smith    TSGetSLES() does not work for integrators that do not use SLES;
1006d763cef2SBarry Smith    in this case TSGetSLES() returns PETSC_NULL in sles.
1007d763cef2SBarry Smith 
1008d763cef2SBarry Smith    Level: beginner
1009d763cef2SBarry Smith 
1010d763cef2SBarry Smith .keywords: timestep, get, SLES
1011d763cef2SBarry Smith @*/
1012d763cef2SBarry Smith int TSGetSLES(TS ts,SLES *sles)
1013d763cef2SBarry Smith {
1014d763cef2SBarry Smith   PetscFunctionBegin;
1015d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
101629bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1017d763cef2SBarry Smith   *sles = ts->sles;
1018d763cef2SBarry Smith   PetscFunctionReturn(0);
1019d763cef2SBarry Smith }
1020d763cef2SBarry Smith 
1021d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1022d763cef2SBarry Smith 
10234a2ae208SSatish Balay #undef __FUNCT__
10244a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1025d763cef2SBarry Smith /*@
1026d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1027d763cef2SBarry Smith    maximum time for iteration.
1028d763cef2SBarry Smith 
1029d763cef2SBarry Smith    Collective on TS
1030d763cef2SBarry Smith 
1031d763cef2SBarry Smith    Input Parameters:
1032d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1033d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1034d763cef2SBarry Smith -  maxtime - final time to iterate to
1035d763cef2SBarry Smith 
1036d763cef2SBarry Smith    Options Database Keys:
1037d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1038d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1039d763cef2SBarry Smith 
1040d763cef2SBarry Smith    Notes:
1041d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1042d763cef2SBarry Smith 
1043d763cef2SBarry Smith    Level: intermediate
1044d763cef2SBarry Smith 
1045d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1046d763cef2SBarry Smith @*/
104787828ca2SBarry Smith int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime)
1048d763cef2SBarry Smith {
1049d763cef2SBarry Smith   PetscFunctionBegin;
1050d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1051d763cef2SBarry Smith   ts->max_steps = maxsteps;
1052d763cef2SBarry Smith   ts->max_time  = maxtime;
1053d763cef2SBarry Smith   PetscFunctionReturn(0);
1054d763cef2SBarry Smith }
1055d763cef2SBarry Smith 
10564a2ae208SSatish Balay #undef __FUNCT__
10574a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1058d763cef2SBarry Smith /*@
1059d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1060d763cef2SBarry Smith    for use by the TS routines.
1061d763cef2SBarry Smith 
1062d763cef2SBarry Smith    Collective on TS and Vec
1063d763cef2SBarry Smith 
1064d763cef2SBarry Smith    Input Parameters:
1065d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1066d763cef2SBarry Smith -  x - the solution vector
1067d763cef2SBarry Smith 
1068d763cef2SBarry Smith    Level: beginner
1069d763cef2SBarry Smith 
1070d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1071d763cef2SBarry Smith @*/
1072d763cef2SBarry Smith int TSSetSolution(TS ts,Vec x)
1073d763cef2SBarry Smith {
1074d763cef2SBarry Smith   PetscFunctionBegin;
1075d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1076d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1077d763cef2SBarry Smith   PetscFunctionReturn(0);
1078d763cef2SBarry Smith }
1079d763cef2SBarry Smith 
1080e74ef692SMatthew Knepley #undef __FUNCT__
1081e74ef692SMatthew Knepley #define __FUNCT__ "TSSetRhsBC"
1082000e7ae3SMatthew Knepley /*@
1083000e7ae3SMatthew Knepley   TSSetRhsBC - Sets the function which applies boundary conditions
1084000e7ae3SMatthew Knepley   to the Rhs of each system.
1085000e7ae3SMatthew Knepley 
1086000e7ae3SMatthew Knepley   Collective on TS
1087000e7ae3SMatthew Knepley 
1088000e7ae3SMatthew Knepley   Input Parameters:
1089000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1090000e7ae3SMatthew Knepley - func - The function
1091000e7ae3SMatthew Knepley 
1092000e7ae3SMatthew Knepley   Calling sequence of func:
1093000e7ae3SMatthew Knepley . func (TS ts, Vec rhs, void *ctx);
1094000e7ae3SMatthew Knepley 
1095000e7ae3SMatthew Knepley + rhs - The current rhs vector
1096000e7ae3SMatthew Knepley - ctx - The user-context
1097000e7ae3SMatthew Knepley 
1098000e7ae3SMatthew Knepley   Level: intermediate
1099000e7ae3SMatthew Knepley 
1100000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1101000e7ae3SMatthew Knepley @*/
1102000e7ae3SMatthew Knepley int TSSetRhsBC(TS ts, int (*func)(TS, Vec, void *))
1103000e7ae3SMatthew Knepley {
1104000e7ae3SMatthew Knepley   PetscFunctionBegin;
1105000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1106000e7ae3SMatthew Knepley   ts->ops->applyrhsbc = func;
1107000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1108000e7ae3SMatthew Knepley }
1109000e7ae3SMatthew Knepley 
1110e74ef692SMatthew Knepley #undef __FUNCT__
1111e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultRhsBC"
1112000e7ae3SMatthew Knepley /*@
1113000e7ae3SMatthew Knepley   TSDefaultRhsBC - The default boundary condition function which does nothing.
1114000e7ae3SMatthew Knepley 
1115000e7ae3SMatthew Knepley   Collective on TS
1116000e7ae3SMatthew Knepley 
1117000e7ae3SMatthew Knepley   Input Parameters:
1118000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1119000e7ae3SMatthew Knepley . rhs - The Rhs
1120000e7ae3SMatthew Knepley - ctx - The user-context
1121000e7ae3SMatthew Knepley 
1122000e7ae3SMatthew Knepley   Level: developer
1123000e7ae3SMatthew Knepley 
1124000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1125000e7ae3SMatthew Knepley @*/
1126000e7ae3SMatthew Knepley int TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1127000e7ae3SMatthew Knepley {
1128000e7ae3SMatthew Knepley   PetscFunctionBegin;
1129000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1130000e7ae3SMatthew Knepley }
1131000e7ae3SMatthew Knepley 
1132e74ef692SMatthew Knepley #undef __FUNCT__
1133e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSystemMatrixBC"
1134000e7ae3SMatthew Knepley /*@
1135000e7ae3SMatthew Knepley   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1136000e7ae3SMatthew Knepley   to the system matrix and preconditioner of each system.
1137000e7ae3SMatthew Knepley 
1138000e7ae3SMatthew Knepley   Collective on TS
1139000e7ae3SMatthew Knepley 
1140000e7ae3SMatthew Knepley   Input Parameters:
1141000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1142000e7ae3SMatthew Knepley - func - The function
1143000e7ae3SMatthew Knepley 
1144000e7ae3SMatthew Knepley   Calling sequence of func:
1145000e7ae3SMatthew Knepley . func (TS ts, Mat A, Mat B, void *ctx);
1146000e7ae3SMatthew Knepley 
1147000e7ae3SMatthew Knepley + A   - The current system matrix
1148000e7ae3SMatthew Knepley . B   - The current preconditioner
1149000e7ae3SMatthew Knepley - ctx - The user-context
1150000e7ae3SMatthew Knepley 
1151000e7ae3SMatthew Knepley   Level: intermediate
1152000e7ae3SMatthew Knepley 
1153000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1154000e7ae3SMatthew Knepley @*/
1155000e7ae3SMatthew Knepley int TSSetSystemMatrixBC(TS ts, int (*func)(TS, Mat, Mat, void *))
1156000e7ae3SMatthew Knepley {
1157000e7ae3SMatthew Knepley   PetscFunctionBegin;
1158000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1159000e7ae3SMatthew Knepley   ts->ops->applymatrixbc = func;
1160000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1161000e7ae3SMatthew Knepley }
1162000e7ae3SMatthew Knepley 
1163e74ef692SMatthew Knepley #undef __FUNCT__
1164e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSystemMatrixBC"
1165000e7ae3SMatthew Knepley /*@
1166000e7ae3SMatthew Knepley   TSDefaultSystemMatrixBC - The default boundary condition function which
1167000e7ae3SMatthew Knepley   does nothing.
1168000e7ae3SMatthew Knepley 
1169000e7ae3SMatthew Knepley   Collective on TS
1170000e7ae3SMatthew Knepley 
1171000e7ae3SMatthew Knepley   Input Parameters:
1172000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1173000e7ae3SMatthew Knepley . A   - The system matrix
1174000e7ae3SMatthew Knepley . B   - The preconditioner
1175000e7ae3SMatthew Knepley - ctx - The user-context
1176000e7ae3SMatthew Knepley 
1177000e7ae3SMatthew Knepley   Level: developer
1178000e7ae3SMatthew Knepley 
1179000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1180000e7ae3SMatthew Knepley @*/
1181000e7ae3SMatthew Knepley int TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1182000e7ae3SMatthew Knepley {
1183000e7ae3SMatthew Knepley   PetscFunctionBegin;
1184000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1185000e7ae3SMatthew Knepley }
1186000e7ae3SMatthew Knepley 
1187e74ef692SMatthew Knepley #undef __FUNCT__
1188e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSolutionBC"
1189000e7ae3SMatthew Knepley /*@
1190000e7ae3SMatthew Knepley   TSSetSolutionBC - Sets the function which applies boundary conditions
1191000e7ae3SMatthew Knepley   to the solution of each system. This is necessary in nonlinear systems
1192000e7ae3SMatthew Knepley   which time dependent boundary conditions.
1193000e7ae3SMatthew Knepley 
1194000e7ae3SMatthew Knepley   Collective on TS
1195000e7ae3SMatthew Knepley 
1196000e7ae3SMatthew Knepley   Input Parameters:
1197000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1198000e7ae3SMatthew Knepley - func - The function
1199000e7ae3SMatthew Knepley 
1200000e7ae3SMatthew Knepley   Calling sequence of func:
1201000e7ae3SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
1202000e7ae3SMatthew Knepley 
1203000e7ae3SMatthew Knepley + sol - The current solution vector
1204000e7ae3SMatthew Knepley - ctx - The user-context
1205000e7ae3SMatthew Knepley 
1206000e7ae3SMatthew Knepley   Level: intermediate
1207000e7ae3SMatthew Knepley 
1208000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1209000e7ae3SMatthew Knepley @*/
1210000e7ae3SMatthew Knepley int TSSetSolutionBC(TS ts, int (*func)(TS, Vec, void *))
1211000e7ae3SMatthew Knepley {
1212000e7ae3SMatthew Knepley   PetscFunctionBegin;
1213000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1214000e7ae3SMatthew Knepley   ts->ops->applysolbc = func;
1215000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1216000e7ae3SMatthew Knepley }
1217000e7ae3SMatthew Knepley 
1218e74ef692SMatthew Knepley #undef __FUNCT__
1219e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSolutionBC"
1220000e7ae3SMatthew Knepley /*@
1221000e7ae3SMatthew Knepley   TSDefaultSolutionBC - The default boundary condition function which
1222000e7ae3SMatthew Knepley   does nothing.
1223000e7ae3SMatthew Knepley 
1224000e7ae3SMatthew Knepley   Collective on TS
1225000e7ae3SMatthew Knepley 
1226000e7ae3SMatthew Knepley   Input Parameters:
1227000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1228000e7ae3SMatthew Knepley . sol - The solution
1229000e7ae3SMatthew Knepley - ctx - The user-context
1230000e7ae3SMatthew Knepley 
1231000e7ae3SMatthew Knepley   Level: developer
1232000e7ae3SMatthew Knepley 
1233000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1234000e7ae3SMatthew Knepley @*/
1235000e7ae3SMatthew Knepley int TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1236000e7ae3SMatthew Knepley {
1237000e7ae3SMatthew Knepley   PetscFunctionBegin;
1238000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1239000e7ae3SMatthew Knepley }
1240000e7ae3SMatthew Knepley 
1241e74ef692SMatthew Knepley #undef __FUNCT__
1242e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1243000e7ae3SMatthew Knepley /*@
1244000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1245000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1246000e7ae3SMatthew Knepley 
1247000e7ae3SMatthew Knepley   Collective on TS
1248000e7ae3SMatthew Knepley 
1249000e7ae3SMatthew Knepley   Input Parameters:
1250000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1251000e7ae3SMatthew Knepley - func - The function
1252000e7ae3SMatthew Knepley 
1253000e7ae3SMatthew Knepley   Calling sequence of func:
1254000e7ae3SMatthew Knepley . func (TS ts);
1255000e7ae3SMatthew Knepley 
1256000e7ae3SMatthew Knepley   Level: intermediate
1257000e7ae3SMatthew Knepley 
1258000e7ae3SMatthew Knepley .keywords: TS, timestep
1259000e7ae3SMatthew Knepley @*/
1260000e7ae3SMatthew Knepley int TSSetPreStep(TS ts, int (*func)(TS))
1261000e7ae3SMatthew Knepley {
1262000e7ae3SMatthew Knepley   PetscFunctionBegin;
1263000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1264000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1265000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1266000e7ae3SMatthew Knepley }
1267000e7ae3SMatthew Knepley 
1268e74ef692SMatthew Knepley #undef __FUNCT__
1269e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1270000e7ae3SMatthew Knepley /*@
1271000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1272000e7ae3SMatthew Knepley 
1273000e7ae3SMatthew Knepley   Collective on TS
1274000e7ae3SMatthew Knepley 
1275000e7ae3SMatthew Knepley   Input Parameters:
1276000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1277000e7ae3SMatthew Knepley 
1278000e7ae3SMatthew Knepley   Level: developer
1279000e7ae3SMatthew Knepley 
1280000e7ae3SMatthew Knepley .keywords: TS, timestep
1281000e7ae3SMatthew Knepley @*/
1282000e7ae3SMatthew Knepley int TSDefaultPreStep(TS ts)
1283000e7ae3SMatthew Knepley {
1284000e7ae3SMatthew Knepley   PetscFunctionBegin;
1285000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1286000e7ae3SMatthew Knepley }
1287000e7ae3SMatthew Knepley 
1288e74ef692SMatthew Knepley #undef __FUNCT__
1289e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1290000e7ae3SMatthew Knepley /*@
1291000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1292000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1293000e7ae3SMatthew Knepley   the time step.
1294000e7ae3SMatthew Knepley 
1295000e7ae3SMatthew Knepley   Collective on TS
1296000e7ae3SMatthew Knepley 
1297000e7ae3SMatthew Knepley   Input Parameters:
1298000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1299000e7ae3SMatthew Knepley - func - The function
1300000e7ae3SMatthew Knepley 
1301000e7ae3SMatthew Knepley   Calling sequence of func:
1302000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1303000e7ae3SMatthew Knepley 
1304000e7ae3SMatthew Knepley + t   - The current time
1305000e7ae3SMatthew Knepley - dt  - The current time step
1306000e7ae3SMatthew Knepley 
1307000e7ae3SMatthew Knepley   Level: intermediate
1308000e7ae3SMatthew Knepley 
1309000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1310000e7ae3SMatthew Knepley @*/
131189cef881SMatthew Knepley int TSSetUpdate(TS ts, int (*func)(TS, PetscReal, PetscReal *))
1312000e7ae3SMatthew Knepley {
1313000e7ae3SMatthew Knepley   PetscFunctionBegin;
1314000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1315000e7ae3SMatthew Knepley   ts->ops->update = func;
1316000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1317000e7ae3SMatthew Knepley }
1318000e7ae3SMatthew Knepley 
1319e74ef692SMatthew Knepley #undef __FUNCT__
1320e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1321000e7ae3SMatthew Knepley /*@
1322000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1323000e7ae3SMatthew Knepley 
1324000e7ae3SMatthew Knepley   Collective on TS
1325000e7ae3SMatthew Knepley 
1326000e7ae3SMatthew Knepley   Input Parameters:
1327000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1328000e7ae3SMatthew Knepley - t   - The current time
1329000e7ae3SMatthew Knepley 
1330000e7ae3SMatthew Knepley   Output Parameters:
1331000e7ae3SMatthew Knepley . dt  - The current time step
1332000e7ae3SMatthew Knepley 
1333000e7ae3SMatthew Knepley   Level: developer
1334000e7ae3SMatthew Knepley 
1335000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1336000e7ae3SMatthew Knepley @*/
133789cef881SMatthew Knepley int TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1338000e7ae3SMatthew Knepley {
1339000e7ae3SMatthew Knepley   PetscFunctionBegin;
1340000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1341000e7ae3SMatthew Knepley }
1342000e7ae3SMatthew Knepley 
1343e74ef692SMatthew Knepley #undef __FUNCT__
1344e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1345000e7ae3SMatthew Knepley /*@
1346000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1347000e7ae3SMatthew Knepley   called once at the end of time stepping.
1348000e7ae3SMatthew Knepley 
1349000e7ae3SMatthew Knepley   Collective on TS
1350000e7ae3SMatthew Knepley 
1351000e7ae3SMatthew Knepley   Input Parameters:
1352000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1353000e7ae3SMatthew Knepley - func - The function
1354000e7ae3SMatthew Knepley 
1355000e7ae3SMatthew Knepley   Calling sequence of func:
1356000e7ae3SMatthew Knepley . func (TS ts);
1357000e7ae3SMatthew Knepley 
1358000e7ae3SMatthew Knepley   Level: intermediate
1359000e7ae3SMatthew Knepley 
1360000e7ae3SMatthew Knepley .keywords: TS, timestep
1361000e7ae3SMatthew Knepley @*/
1362000e7ae3SMatthew Knepley int TSSetPostStep(TS ts, int (*func)(TS))
1363000e7ae3SMatthew Knepley {
1364000e7ae3SMatthew Knepley   PetscFunctionBegin;
1365000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1366000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1367000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1368000e7ae3SMatthew Knepley }
1369000e7ae3SMatthew Knepley 
1370e74ef692SMatthew Knepley #undef __FUNCT__
1371e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1372000e7ae3SMatthew Knepley /*@
1373000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1374000e7ae3SMatthew Knepley 
1375000e7ae3SMatthew Knepley   Collective on TS
1376000e7ae3SMatthew Knepley 
1377000e7ae3SMatthew Knepley   Input Parameters:
1378000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1379000e7ae3SMatthew Knepley 
1380000e7ae3SMatthew Knepley   Level: developer
1381000e7ae3SMatthew Knepley 
1382000e7ae3SMatthew Knepley .keywords: TS, timestep
1383000e7ae3SMatthew Knepley @*/
1384000e7ae3SMatthew Knepley int TSDefaultPostStep(TS ts)
1385000e7ae3SMatthew Knepley {
1386000e7ae3SMatthew Knepley   PetscFunctionBegin;
1387000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1388000e7ae3SMatthew Knepley }
1389000e7ae3SMatthew Knepley 
1390d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1391d763cef2SBarry Smith 
13924a2ae208SSatish Balay #undef __FUNCT__
13934a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
1394d763cef2SBarry Smith /*@C
1395d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1396d763cef2SBarry Smith    timestep to display the iteration's  progress.
1397d763cef2SBarry Smith 
1398d763cef2SBarry Smith    Collective on TS
1399d763cef2SBarry Smith 
1400d763cef2SBarry Smith    Input Parameters:
1401d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1402d763cef2SBarry Smith .  func - monitoring routine
1403329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1404b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1405b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1406b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1407d763cef2SBarry Smith 
1408d763cef2SBarry Smith    Calling sequence of func:
140987828ca2SBarry Smith $    int func(TS ts,int steps,PetscReal time,Vec x,void *mctx)
1410d763cef2SBarry Smith 
1411d763cef2SBarry Smith +    ts - the TS context
1412d763cef2SBarry Smith .    steps - iteration number
14131f06c33eSBarry Smith .    time - current time
1414d763cef2SBarry Smith .    x - current iterate
1415d763cef2SBarry Smith -    mctx - [optional] monitoring context
1416d763cef2SBarry Smith 
1417d763cef2SBarry Smith    Notes:
1418d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1419d763cef2SBarry Smith    already has been loaded.
1420d763cef2SBarry Smith 
1421d763cef2SBarry Smith    Level: intermediate
1422d763cef2SBarry Smith 
1423d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1424d763cef2SBarry Smith 
1425d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
1426d763cef2SBarry Smith @*/
142787828ca2SBarry Smith int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*))
1428d763cef2SBarry Smith {
1429d763cef2SBarry Smith   PetscFunctionBegin;
1430d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1431d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
143229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1433d763cef2SBarry Smith   }
1434d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1435329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1436d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1437d763cef2SBarry Smith   PetscFunctionReturn(0);
1438d763cef2SBarry Smith }
1439d763cef2SBarry Smith 
14404a2ae208SSatish Balay #undef __FUNCT__
14414a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
1442d763cef2SBarry Smith /*@C
1443d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1444d763cef2SBarry Smith 
1445d763cef2SBarry Smith    Collective on TS
1446d763cef2SBarry Smith 
1447d763cef2SBarry Smith    Input Parameters:
1448d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1449d763cef2SBarry Smith 
1450d763cef2SBarry Smith    Notes:
1451d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1452d763cef2SBarry Smith 
1453d763cef2SBarry Smith    Level: intermediate
1454d763cef2SBarry Smith 
1455d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1456d763cef2SBarry Smith 
1457d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
1458d763cef2SBarry Smith @*/
1459d763cef2SBarry Smith int TSClearMonitor(TS ts)
1460d763cef2SBarry Smith {
1461d763cef2SBarry Smith   PetscFunctionBegin;
1462d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1463d763cef2SBarry Smith   ts->numbermonitors = 0;
1464d763cef2SBarry Smith   PetscFunctionReturn(0);
1465d763cef2SBarry Smith }
1466d763cef2SBarry Smith 
14674a2ae208SSatish Balay #undef __FUNCT__
14684a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
146987828ca2SBarry Smith int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
1470d763cef2SBarry Smith {
1471d132466eSBarry Smith   int ierr;
1472d132466eSBarry Smith 
1473d763cef2SBarry Smith   PetscFunctionBegin;
1474142b95e3SSatish Balay   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1475d763cef2SBarry Smith   PetscFunctionReturn(0);
1476d763cef2SBarry Smith }
1477d763cef2SBarry Smith 
14784a2ae208SSatish Balay #undef __FUNCT__
14794a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1480d763cef2SBarry Smith /*@
1481d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1482d763cef2SBarry Smith 
1483d763cef2SBarry Smith    Collective on TS
1484d763cef2SBarry Smith 
1485d763cef2SBarry Smith    Input Parameter:
1486d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1487d763cef2SBarry Smith 
1488d763cef2SBarry Smith    Output Parameters:
1489d763cef2SBarry Smith +  steps - number of iterations until termination
1490142b95e3SSatish Balay -  ptime - time until termination
1491d763cef2SBarry Smith 
1492d763cef2SBarry Smith    Level: beginner
1493d763cef2SBarry Smith 
1494d763cef2SBarry Smith .keywords: TS, timestep, solve
1495d763cef2SBarry Smith 
1496d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1497d763cef2SBarry Smith @*/
149887828ca2SBarry Smith int TSStep(TS ts,int *steps,PetscReal *ptime)
1499d763cef2SBarry Smith {
1500d405a339SMatthew Knepley   PetscTruth opt;
1501f1af5d2fSBarry Smith   int        ierr;
1502d763cef2SBarry Smith 
1503d763cef2SBarry Smith   PetscFunctionBegin;
1504d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1505d405a339SMatthew Knepley   if (!ts->setupcalled) {
1506d405a339SMatthew Knepley     ierr = TSSetUp(ts);                                                                                   CHKERRQ(ierr);
1507d405a339SMatthew Knepley   }
1508d405a339SMatthew Knepley 
1509d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);                                                        CHKERRQ(ierr);
1510d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);                                                                         CHKERRQ(ierr);
1511000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);                                                              CHKERRQ(ierr);
1512d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);                                                                        CHKERRQ(ierr);
1513d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);                                                          CHKERRQ(ierr);
1514d405a339SMatthew Knepley 
1515*4bb05414SBarry Smith   if (!PetscPreLoadingOn) {
1516*4bb05414SBarry Smith     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1517d405a339SMatthew Knepley   }
1518d763cef2SBarry Smith   PetscFunctionReturn(0);
1519d763cef2SBarry Smith }
1520d763cef2SBarry Smith 
15214a2ae208SSatish Balay #undef __FUNCT__
15224a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1523d763cef2SBarry Smith /*
1524d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1525d763cef2SBarry Smith */
152687828ca2SBarry Smith int TSMonitor(TS ts,int step,PetscReal ptime,Vec x)
1527d763cef2SBarry Smith {
1528d763cef2SBarry Smith   int i,ierr,n = ts->numbermonitors;
1529d763cef2SBarry Smith 
1530d763cef2SBarry Smith   PetscFunctionBegin;
1531d763cef2SBarry Smith   for (i=0; i<n; i++) {
1532142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1533d763cef2SBarry Smith   }
1534d763cef2SBarry Smith   PetscFunctionReturn(0);
1535d763cef2SBarry Smith }
1536d763cef2SBarry Smith 
1537d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1538d763cef2SBarry Smith 
15394a2ae208SSatish Balay #undef __FUNCT__
15404a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1541d763cef2SBarry Smith /*@C
1542d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1543d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1544d763cef2SBarry Smith 
1545d763cef2SBarry Smith    Collective on TS
1546d763cef2SBarry Smith 
1547d763cef2SBarry Smith    Input Parameters:
1548d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1549d763cef2SBarry Smith .  label - the title to put in the title bar
15507c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1551d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1552d763cef2SBarry Smith 
1553d763cef2SBarry Smith    Output Parameter:
1554d763cef2SBarry Smith .  draw - the drawing context
1555d763cef2SBarry Smith 
1556d763cef2SBarry Smith    Options Database Key:
1557d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1558d763cef2SBarry Smith 
1559d763cef2SBarry Smith    Notes:
1560b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1561d763cef2SBarry Smith 
1562d763cef2SBarry Smith    Level: intermediate
1563d763cef2SBarry Smith 
15647c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1565d763cef2SBarry Smith 
1566d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
15677c922b88SBarry Smith 
1568d763cef2SBarry Smith @*/
1569b0a32e0cSBarry Smith int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1570d763cef2SBarry Smith {
1571b0a32e0cSBarry Smith   PetscDraw win;
1572d763cef2SBarry Smith   int       ierr;
1573d763cef2SBarry Smith 
1574d763cef2SBarry Smith   PetscFunctionBegin;
1575b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1576b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1577b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1578b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1579d763cef2SBarry Smith 
1580b0a32e0cSBarry Smith   PetscLogObjectParent(*draw,win);
1581d763cef2SBarry Smith   PetscFunctionReturn(0);
1582d763cef2SBarry Smith }
1583d763cef2SBarry Smith 
15844a2ae208SSatish Balay #undef __FUNCT__
15854a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
158687828ca2SBarry Smith int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx)
1587d763cef2SBarry Smith {
1588b0a32e0cSBarry Smith   PetscDrawLG lg = (PetscDrawLG) monctx;
158987828ca2SBarry Smith   PetscReal      x,y = ptime;
1590d763cef2SBarry Smith   int         ierr;
1591d763cef2SBarry Smith 
1592d763cef2SBarry Smith   PetscFunctionBegin;
15937c922b88SBarry Smith   if (!monctx) {
15947c922b88SBarry Smith     MPI_Comm    comm;
1595b0a32e0cSBarry Smith     PetscViewer viewer;
15967c922b88SBarry Smith 
15977c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1598b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1599b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
16007c922b88SBarry Smith   }
16017c922b88SBarry Smith 
1602b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
160387828ca2SBarry Smith   x = (PetscReal)n;
1604b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1605d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1606b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1607d763cef2SBarry Smith   }
1608d763cef2SBarry Smith   PetscFunctionReturn(0);
1609d763cef2SBarry Smith }
1610d763cef2SBarry Smith 
16114a2ae208SSatish Balay #undef __FUNCT__
16124a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1613d763cef2SBarry Smith /*@C
1614d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1615d763cef2SBarry Smith    with TSLGMonitorCreate().
1616d763cef2SBarry Smith 
1617b0a32e0cSBarry Smith    Collective on PetscDrawLG
1618d763cef2SBarry Smith 
1619d763cef2SBarry Smith    Input Parameter:
1620d763cef2SBarry Smith .  draw - the drawing context
1621d763cef2SBarry Smith 
1622d763cef2SBarry Smith    Level: intermediate
1623d763cef2SBarry Smith 
1624d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1625d763cef2SBarry Smith 
1626d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1627d763cef2SBarry Smith @*/
1628b0a32e0cSBarry Smith int TSLGMonitorDestroy(PetscDrawLG drawlg)
1629d763cef2SBarry Smith {
1630b0a32e0cSBarry Smith   PetscDraw draw;
1631d763cef2SBarry Smith   int       ierr;
1632d763cef2SBarry Smith 
1633d763cef2SBarry Smith   PetscFunctionBegin;
1634b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1635b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1636b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1637d763cef2SBarry Smith   PetscFunctionReturn(0);
1638d763cef2SBarry Smith }
1639d763cef2SBarry Smith 
16404a2ae208SSatish Balay #undef __FUNCT__
16414a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1642d763cef2SBarry Smith /*@
1643d763cef2SBarry Smith    TSGetTime - Gets the current time.
1644d763cef2SBarry Smith 
1645d763cef2SBarry Smith    Not Collective
1646d763cef2SBarry Smith 
1647d763cef2SBarry Smith    Input Parameter:
1648d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1649d763cef2SBarry Smith 
1650d763cef2SBarry Smith    Output Parameter:
1651d763cef2SBarry Smith .  t  - the current time
1652d763cef2SBarry Smith 
1653d763cef2SBarry Smith    Contributed by: Matthew Knepley
1654d763cef2SBarry Smith 
1655d763cef2SBarry Smith    Level: beginner
1656d763cef2SBarry Smith 
1657d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1658d763cef2SBarry Smith 
1659d763cef2SBarry Smith .keywords: TS, get, time
1660d763cef2SBarry Smith @*/
166187828ca2SBarry Smith int TSGetTime(TS ts,PetscReal* t)
1662d763cef2SBarry Smith {
1663d763cef2SBarry Smith   PetscFunctionBegin;
1664d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1665d763cef2SBarry Smith   *t = ts->ptime;
1666d763cef2SBarry Smith   PetscFunctionReturn(0);
1667d763cef2SBarry Smith }
1668d763cef2SBarry Smith 
16694a2ae208SSatish Balay #undef __FUNCT__
16704a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1671d763cef2SBarry Smith /*@C
1672d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1673d763cef2SBarry Smith    TS options in the database.
1674d763cef2SBarry Smith 
1675d763cef2SBarry Smith    Collective on TS
1676d763cef2SBarry Smith 
1677d763cef2SBarry Smith    Input Parameter:
1678d763cef2SBarry Smith +  ts     - The TS context
1679d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1680d763cef2SBarry Smith 
1681d763cef2SBarry Smith    Notes:
1682d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1683d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1684d763cef2SBarry Smith    hyphen.
1685d763cef2SBarry Smith 
1686d763cef2SBarry Smith    Contributed by: Matthew Knepley
1687d763cef2SBarry Smith 
1688d763cef2SBarry Smith    Level: advanced
1689d763cef2SBarry Smith 
1690d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1691d763cef2SBarry Smith 
1692d763cef2SBarry Smith .seealso: TSSetFromOptions()
1693d763cef2SBarry Smith 
1694d763cef2SBarry Smith @*/
1695d763cef2SBarry Smith int TSSetOptionsPrefix(TS ts,char *prefix)
1696d763cef2SBarry Smith {
1697d763cef2SBarry Smith   int ierr;
1698d763cef2SBarry Smith 
1699d763cef2SBarry Smith   PetscFunctionBegin;
1700d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1701d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1702d763cef2SBarry Smith   switch(ts->problem_type) {
1703d763cef2SBarry Smith     case TS_NONLINEAR:
1704d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1705d763cef2SBarry Smith       break;
1706d763cef2SBarry Smith     case TS_LINEAR:
1707d763cef2SBarry Smith       ierr = SLESSetOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr);
1708d763cef2SBarry Smith       break;
1709d763cef2SBarry Smith   }
1710d763cef2SBarry Smith   PetscFunctionReturn(0);
1711d763cef2SBarry Smith }
1712d763cef2SBarry Smith 
1713d763cef2SBarry Smith 
17144a2ae208SSatish Balay #undef __FUNCT__
17154a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1716d763cef2SBarry Smith /*@C
1717d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1718d763cef2SBarry Smith    TS options in the database.
1719d763cef2SBarry Smith 
1720d763cef2SBarry Smith    Collective on TS
1721d763cef2SBarry Smith 
1722d763cef2SBarry Smith    Input Parameter:
1723d763cef2SBarry Smith +  ts     - The TS context
1724d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1725d763cef2SBarry Smith 
1726d763cef2SBarry Smith    Notes:
1727d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1728d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1729d763cef2SBarry Smith    hyphen.
1730d763cef2SBarry Smith 
1731d763cef2SBarry Smith    Contributed by: Matthew Knepley
1732d763cef2SBarry Smith 
1733d763cef2SBarry Smith    Level: advanced
1734d763cef2SBarry Smith 
1735d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1736d763cef2SBarry Smith 
1737d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1738d763cef2SBarry Smith 
1739d763cef2SBarry Smith @*/
1740d763cef2SBarry Smith int TSAppendOptionsPrefix(TS ts,char *prefix)
1741d763cef2SBarry Smith {
1742d763cef2SBarry Smith   int ierr;
1743d763cef2SBarry Smith 
1744d763cef2SBarry Smith   PetscFunctionBegin;
1745d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1746d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1747d763cef2SBarry Smith   switch(ts->problem_type) {
1748d763cef2SBarry Smith     case TS_NONLINEAR:
1749d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1750d763cef2SBarry Smith       break;
1751d763cef2SBarry Smith     case TS_LINEAR:
1752d763cef2SBarry Smith       ierr = SLESAppendOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr);
1753d763cef2SBarry Smith       break;
1754d763cef2SBarry Smith   }
1755d763cef2SBarry Smith   PetscFunctionReturn(0);
1756d763cef2SBarry Smith }
1757d763cef2SBarry Smith 
17584a2ae208SSatish Balay #undef __FUNCT__
17594a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1760d763cef2SBarry Smith /*@C
1761d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1762d763cef2SBarry Smith    TS options in the database.
1763d763cef2SBarry Smith 
1764d763cef2SBarry Smith    Not Collective
1765d763cef2SBarry Smith 
1766d763cef2SBarry Smith    Input Parameter:
1767d763cef2SBarry Smith .  ts - The TS context
1768d763cef2SBarry Smith 
1769d763cef2SBarry Smith    Output Parameter:
1770d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1771d763cef2SBarry Smith 
1772d763cef2SBarry Smith    Contributed by: Matthew Knepley
1773d763cef2SBarry Smith 
1774d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1775d763cef2SBarry Smith    sufficient length to hold the prefix.
1776d763cef2SBarry Smith 
1777d763cef2SBarry Smith    Level: intermediate
1778d763cef2SBarry Smith 
1779d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1780d763cef2SBarry Smith 
1781d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1782d763cef2SBarry Smith @*/
1783d763cef2SBarry Smith int TSGetOptionsPrefix(TS ts,char **prefix)
1784d763cef2SBarry Smith {
1785d763cef2SBarry Smith   int ierr;
1786d763cef2SBarry Smith 
1787d763cef2SBarry Smith   PetscFunctionBegin;
1788d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1789d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1790d763cef2SBarry Smith   PetscFunctionReturn(0);
1791d763cef2SBarry Smith }
1792d763cef2SBarry Smith 
17934a2ae208SSatish Balay #undef __FUNCT__
17944a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1795d763cef2SBarry Smith /*@C
1796d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1797d763cef2SBarry Smith 
1798d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1799d763cef2SBarry Smith 
1800d763cef2SBarry Smith    Input Parameter:
1801d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1802d763cef2SBarry Smith 
1803d763cef2SBarry Smith    Output Parameters:
1804d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1805d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1806d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1807d763cef2SBarry Smith 
1808d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1809d763cef2SBarry Smith 
1810d763cef2SBarry Smith    Contributed by: Matthew Knepley
1811d763cef2SBarry Smith 
1812d763cef2SBarry Smith    Level: intermediate
1813d763cef2SBarry Smith 
1814d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1815d763cef2SBarry Smith 
1816d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1817d763cef2SBarry Smith 
1818d763cef2SBarry Smith @*/
1819d763cef2SBarry Smith int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1820d763cef2SBarry Smith {
1821d763cef2SBarry Smith   PetscFunctionBegin;
1822d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1823d763cef2SBarry Smith   if (A)   *A = ts->A;
1824d763cef2SBarry Smith   if (M)   *M = ts->B;
1825d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1826d763cef2SBarry Smith   PetscFunctionReturn(0);
1827d763cef2SBarry Smith }
1828d763cef2SBarry Smith 
18294a2ae208SSatish Balay #undef __FUNCT__
18304a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1831d763cef2SBarry Smith /*@C
1832d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1833d763cef2SBarry Smith 
1834d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1835d763cef2SBarry Smith 
1836d763cef2SBarry Smith    Input Parameter:
1837d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1838d763cef2SBarry Smith 
1839d763cef2SBarry Smith    Output Parameters:
1840d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1841d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1842d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1843d763cef2SBarry Smith 
1844d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1845d763cef2SBarry Smith 
1846d763cef2SBarry Smith    Contributed by: Matthew Knepley
1847d763cef2SBarry Smith 
1848d763cef2SBarry Smith    Level: intermediate
1849d763cef2SBarry Smith 
1850d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1851d763cef2SBarry Smith 
1852d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1853d763cef2SBarry Smith @*/
1854d763cef2SBarry Smith int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1855d763cef2SBarry Smith {
1856d763cef2SBarry Smith   int ierr;
1857d763cef2SBarry Smith 
1858d763cef2SBarry Smith   PetscFunctionBegin;
1859d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1860d763cef2SBarry Smith   PetscFunctionReturn(0);
1861d763cef2SBarry Smith }
1862d763cef2SBarry Smith 
18631713a123SBarry Smith #undef __FUNCT__
18641713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
18651713a123SBarry Smith /*@C
18661713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
18671713a123SBarry Smith    VecView() for the solution at each timestep
18681713a123SBarry Smith 
18691713a123SBarry Smith    Collective on TS
18701713a123SBarry Smith 
18711713a123SBarry Smith    Input Parameters:
18721713a123SBarry Smith +  ts - the TS context
18731713a123SBarry Smith .  step - current time-step
1874142b95e3SSatish Balay .  ptime - current time
18751713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
18761713a123SBarry Smith 
18771713a123SBarry Smith    Level: intermediate
18781713a123SBarry Smith 
18791713a123SBarry Smith .keywords: TS,  vector, monitor, view
18801713a123SBarry Smith 
18811713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
18821713a123SBarry Smith @*/
1883142b95e3SSatish Balay int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
18841713a123SBarry Smith {
18851713a123SBarry Smith   int         ierr;
18861713a123SBarry Smith   PetscViewer viewer = (PetscViewer) dummy;
18871713a123SBarry Smith 
18881713a123SBarry Smith   PetscFunctionBegin;
18891713a123SBarry Smith   if (!viewer) {
18901713a123SBarry Smith     MPI_Comm comm;
18911713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
18921713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
18931713a123SBarry Smith   }
18941713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
18951713a123SBarry Smith   PetscFunctionReturn(0);
18961713a123SBarry Smith }
18971713a123SBarry Smith 
18981713a123SBarry Smith 
18991713a123SBarry Smith 
1900