xref: /petsc/src/ts/interface/ts.c (revision 2fc52814d27bf1f4e71021c1c3ebb532b583ed60)
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 */
543c77886SBarry Smith int TS_COOKIE = 0;
643c77886SBarry Smith int TS_Step = 0, TS_PseudoComputeTimeStep = 0, TS_FunctionEval = 0, TS_JacobianEval = 0;
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;
26*2fc52814SBarry Smith   const 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:
11994b7f48cSBarry Smith     if (ts->ksp != PETSC_NULL) {
12094b7f48cSBarry Smith       ierr = KSPSetFromOptions(ts->ksp);                                                                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 @*/
1501836bdbcSSatish Balay int TSViewFromOptions(TS ts,const char title[])
151bdad233fSMatthew Knepley {
152bdad233fSMatthew Knepley   PetscViewer viewer;
153bdad233fSMatthew Knepley   PetscDraw   draw;
154bdad233fSMatthew Knepley   PetscTruth  opt;
155bdad233fSMatthew Knepley   char        typeName[1024];
156e10c49a3SBarry Smith   char        fileName[PETSC_MAX_PATH_LEN];
157bdad233fSMatthew Knepley   int         len;
158bdad233fSMatthew Knepley   int         ierr;
159bdad233fSMatthew Knepley 
160bdad233fSMatthew Knepley   PetscFunctionBegin;
161bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);                                               CHKERRQ(ierr);
162bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
163bdad233fSMatthew Knepley     ierr = PetscOptionsGetString(ts->prefix, "-ts_view", typeName, 1024, &opt);                           CHKERRQ(ierr);
164bdad233fSMatthew Knepley     ierr = PetscStrlen(typeName, &len);                                                                   CHKERRQ(ierr);
165bdad233fSMatthew Knepley     if (len > 0) {
166bdad233fSMatthew Knepley       ierr = PetscViewerCreate(ts->comm, &viewer);                                                        CHKERRQ(ierr);
167bdad233fSMatthew Knepley       ierr = PetscViewerSetType(viewer, typeName);                                                        CHKERRQ(ierr);
168bdad233fSMatthew Knepley       ierr = PetscOptionsGetString(ts->prefix, "-ts_view_file", fileName, 1024, &opt);                    CHKERRQ(ierr);
169bdad233fSMatthew Knepley       if (opt == PETSC_TRUE) {
170bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, fileName);                                                  CHKERRQ(ierr);
171bdad233fSMatthew Knepley       } else {
172bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, ts->name);                                                  CHKERRQ(ierr);
173bdad233fSMatthew Knepley       }
174bdad233fSMatthew Knepley       ierr = TSView(ts, viewer);                                                                          CHKERRQ(ierr);
175bdad233fSMatthew Knepley       ierr = PetscViewerFlush(viewer);                                                                    CHKERRQ(ierr);
176bdad233fSMatthew Knepley       ierr = PetscViewerDestroy(viewer);                                                                  CHKERRQ(ierr);
177bdad233fSMatthew Knepley     } else {
178bdad233fSMatthew Knepley       ierr = TSView(ts, PETSC_NULL);                                                                      CHKERRQ(ierr);
179bdad233fSMatthew Knepley     }
180bdad233fSMatthew Knepley   }
181bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);                                          CHKERRQ(ierr);
182bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
183bdad233fSMatthew Knepley     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);                                  CHKERRQ(ierr);
184bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);                                                      CHKERRQ(ierr);
185bdad233fSMatthew Knepley     if (title != PETSC_NULL) {
1861836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, (char *)title);                                                      CHKERRQ(ierr);
187bdad233fSMatthew Knepley     } else {
188bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject) ts);                                                           CHKERRQ(ierr) ;
1891836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, ts->name);                                                           CHKERRQ(ierr);
190bdad233fSMatthew Knepley     }
191bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);                                                                            CHKERRQ(ierr);
192bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);                                                                      CHKERRQ(ierr);
193bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);                                                                          CHKERRQ(ierr);
194bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);                                                                    CHKERRQ(ierr);
195bdad233fSMatthew Knepley   }
196bdad233fSMatthew Knepley   PetscFunctionReturn(0);
197bdad233fSMatthew Knepley }
198bdad233fSMatthew Knepley 
199bdad233fSMatthew Knepley #undef __FUNCT__
2004a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
201a7a1495cSBarry Smith /*@
2028c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
203a7a1495cSBarry Smith       set with TSSetRHSJacobian().
204a7a1495cSBarry Smith 
205a7a1495cSBarry Smith    Collective on TS and Vec
206a7a1495cSBarry Smith 
207a7a1495cSBarry Smith    Input Parameters:
208a7a1495cSBarry Smith +  ts - the SNES context
209a7a1495cSBarry Smith .  t - current timestep
210a7a1495cSBarry Smith -  x - input vector
211a7a1495cSBarry Smith 
212a7a1495cSBarry Smith    Output Parameters:
213a7a1495cSBarry Smith +  A - Jacobian matrix
214a7a1495cSBarry Smith .  B - optional preconditioning matrix
215a7a1495cSBarry Smith -  flag - flag indicating matrix structure
216a7a1495cSBarry Smith 
217a7a1495cSBarry Smith    Notes:
218a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
219a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
220a7a1495cSBarry Smith 
22194b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
222a7a1495cSBarry Smith    flag parameter.
223a7a1495cSBarry Smith 
224a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
225a7a1495cSBarry Smith 
226a7a1495cSBarry Smith    Level: developer
227a7a1495cSBarry Smith 
228a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
229a7a1495cSBarry Smith 
23094b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
231a7a1495cSBarry Smith @*/
23287828ca2SBarry Smith int TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
233a7a1495cSBarry Smith {
234a7a1495cSBarry Smith   int ierr;
235a7a1495cSBarry Smith 
236a7a1495cSBarry Smith   PetscFunctionBegin;
237a7a1495cSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
238a7a1495cSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE);
239a7a1495cSBarry Smith   PetscCheckSameComm(ts,X);
240a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
24129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
242a7a1495cSBarry Smith   }
243000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
244d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
245a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
246a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
247000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
248a7a1495cSBarry Smith     PetscStackPop;
249d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
250a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
251a7a1495cSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE);
252a7a1495cSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE);
253ef66eb69SBarry Smith   } else {
254ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256ef66eb69SBarry Smith     if (*A != *B) {
257ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259ef66eb69SBarry Smith     }
260ef66eb69SBarry Smith   }
261a7a1495cSBarry Smith   PetscFunctionReturn(0);
262a7a1495cSBarry Smith }
263a7a1495cSBarry Smith 
2644a2ae208SSatish Balay #undef __FUNCT__
2654a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
266d763cef2SBarry Smith /*
267d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
268d763cef2SBarry Smith 
269d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
270d763cef2SBarry Smith    this routine applies the matrix.
271d763cef2SBarry Smith */
27287828ca2SBarry Smith int TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
273d763cef2SBarry Smith {
274d763cef2SBarry Smith   int ierr;
275d763cef2SBarry Smith 
276d763cef2SBarry Smith   PetscFunctionBegin;
277d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
2787c922b88SBarry Smith   PetscValidHeader(x);
2797c922b88SBarry Smith   PetscValidHeader(y);
280d763cef2SBarry Smith 
281d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
282000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
283d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
284000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
285d763cef2SBarry Smith     PetscStackPop;
2867c922b88SBarry Smith   } else {
287000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
288d763cef2SBarry Smith       MatStructure flg;
289d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
290000e7ae3SMatthew Knepley       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
291d763cef2SBarry Smith       PetscStackPop;
292d763cef2SBarry Smith     }
293d763cef2SBarry Smith     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
2947c922b88SBarry Smith   }
295d763cef2SBarry Smith 
296d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
297d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
298d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
299d763cef2SBarry Smith 
300d763cef2SBarry Smith   PetscFunctionReturn(0);
301d763cef2SBarry Smith }
302d763cef2SBarry Smith 
3034a2ae208SSatish Balay #undef __FUNCT__
3044a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
305d763cef2SBarry Smith /*@C
306d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
307d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
308d763cef2SBarry Smith 
309d763cef2SBarry Smith     Collective on TS
310d763cef2SBarry Smith 
311d763cef2SBarry Smith     Input Parameters:
312d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
313d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
314d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
315d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
316d763cef2SBarry Smith 
317d763cef2SBarry Smith     Calling sequence of func:
31887828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
319d763cef2SBarry Smith 
320d763cef2SBarry Smith +   t - current timestep
321d763cef2SBarry Smith .   u - input vector
322d763cef2SBarry Smith .   F - function vector
323d763cef2SBarry Smith -   ctx - [optional] user-defined function context
324d763cef2SBarry Smith 
325d763cef2SBarry Smith     Important:
326d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
327d763cef2SBarry Smith 
328d763cef2SBarry Smith     Level: beginner
329d763cef2SBarry Smith 
330d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
331d763cef2SBarry Smith 
332d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
333d763cef2SBarry Smith @*/
33487828ca2SBarry Smith int TSSetRHSFunction(TS ts,int (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
335d763cef2SBarry Smith {
336d763cef2SBarry Smith   PetscFunctionBegin;
337d763cef2SBarry Smith 
338d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
339d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
34029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
341d763cef2SBarry Smith   }
342000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
343d763cef2SBarry Smith   ts->funP             = ctx;
344d763cef2SBarry Smith   PetscFunctionReturn(0);
345d763cef2SBarry Smith }
346d763cef2SBarry Smith 
3474a2ae208SSatish Balay #undef __FUNCT__
3484a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
349d763cef2SBarry Smith /*@C
350d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
351d763cef2SBarry Smith    Also sets the location to store A.
352d763cef2SBarry Smith 
353d763cef2SBarry Smith    Collective on TS
354d763cef2SBarry Smith 
355d763cef2SBarry Smith    Input Parameters:
356d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
357d763cef2SBarry Smith .  A   - matrix
358d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
359d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
360d763cef2SBarry Smith          if A is not a function of t.
361d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
362d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
363d763cef2SBarry Smith 
364d763cef2SBarry Smith    Calling sequence of func:
36587828ca2SBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,int *flag,void *ctx);
366d763cef2SBarry Smith 
367d763cef2SBarry Smith +  t - current timestep
368d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
369d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
370d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
37194b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
372d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
373d763cef2SBarry Smith 
374d763cef2SBarry Smith    Notes:
37594b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
376d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
377d763cef2SBarry Smith 
378d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
379d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
380d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
381d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
382d763cef2SBarry Smith    throughout the global iterations.
383d763cef2SBarry Smith 
384d763cef2SBarry Smith    Important:
385d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
386d763cef2SBarry Smith 
387d763cef2SBarry Smith    Level: beginner
388d763cef2SBarry Smith 
389d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
390d763cef2SBarry Smith 
391d763cef2SBarry Smith .seealso: TSSetRHSFunction()
392d763cef2SBarry Smith @*/
39387828ca2SBarry Smith int TSSetRHSMatrix(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
394d763cef2SBarry Smith {
395d763cef2SBarry Smith   PetscFunctionBegin;
396d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
397184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
398184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
399184914b5SBarry Smith   PetscCheckSameComm(ts,A);
400184914b5SBarry Smith   PetscCheckSameComm(ts,B);
401d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
40229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
403d763cef2SBarry Smith   }
404d763cef2SBarry Smith 
405000e7ae3SMatthew Knepley   ts->ops->rhsmatrix = f;
406d763cef2SBarry Smith   ts->jacP           = ctx;
407d763cef2SBarry Smith   ts->A              = A;
408d763cef2SBarry Smith   ts->B              = B;
409d763cef2SBarry Smith 
410d763cef2SBarry Smith   PetscFunctionReturn(0);
411d763cef2SBarry Smith }
412d763cef2SBarry Smith 
4134a2ae208SSatish Balay #undef __FUNCT__
4144a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
415d763cef2SBarry Smith /*@C
416d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
417d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
418d763cef2SBarry Smith 
419d763cef2SBarry Smith    Collective on TS
420d763cef2SBarry Smith 
421d763cef2SBarry Smith    Input Parameters:
422d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
423d763cef2SBarry Smith .  A   - Jacobian matrix
424d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
425d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
426d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
427d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
428d763cef2SBarry Smith 
429d763cef2SBarry Smith    Calling sequence of func:
43087828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
431d763cef2SBarry Smith 
432d763cef2SBarry Smith +  t - current timestep
433d763cef2SBarry Smith .  u - input vector
434d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
435d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
436d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
43794b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
438d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
439d763cef2SBarry Smith 
440d763cef2SBarry Smith    Notes:
44194b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
442d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
443d763cef2SBarry Smith 
444d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
445d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
446d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
447d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
448d763cef2SBarry Smith    throughout the global iterations.
449d763cef2SBarry Smith 
450d763cef2SBarry Smith    Level: beginner
451d763cef2SBarry Smith 
452d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
453d763cef2SBarry Smith 
454d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
455d763cef2SBarry Smith           SNESDefaultComputeJacobianColor()
456d763cef2SBarry Smith 
457d763cef2SBarry Smith @*/
45887828ca2SBarry Smith int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
459d763cef2SBarry Smith {
460d763cef2SBarry Smith   PetscFunctionBegin;
461d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
462184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
463184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
464184914b5SBarry Smith   PetscCheckSameComm(ts,A);
465184914b5SBarry Smith   PetscCheckSameComm(ts,B);
466d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
46729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
468d763cef2SBarry Smith   }
469d763cef2SBarry Smith 
470000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
471d763cef2SBarry Smith   ts->jacP             = ctx;
472d763cef2SBarry Smith   ts->A                = A;
473d763cef2SBarry Smith   ts->B                = B;
474d763cef2SBarry Smith   PetscFunctionReturn(0);
475d763cef2SBarry Smith }
476d763cef2SBarry Smith 
4774a2ae208SSatish Balay #undef __FUNCT__
4784a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
479d763cef2SBarry Smith /*
480d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
481d763cef2SBarry Smith 
482d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
483d763cef2SBarry Smith    this routine applies the matrix.
484d763cef2SBarry Smith */
48587828ca2SBarry Smith int TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
486d763cef2SBarry Smith {
487d763cef2SBarry Smith   int ierr;
488d763cef2SBarry Smith 
489d763cef2SBarry Smith   PetscFunctionBegin;
490d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
491d763cef2SBarry Smith   PetscValidHeader(x);
492184914b5SBarry Smith   PetscCheckSameComm(ts,x);
493d763cef2SBarry Smith 
494000e7ae3SMatthew Knepley   if (ts->ops->rhsbc) {
495d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
496000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
497d763cef2SBarry Smith     PetscStackPop;
498d763cef2SBarry Smith     PetscFunctionReturn(0);
499d763cef2SBarry Smith   }
500d763cef2SBarry Smith 
501d763cef2SBarry Smith   PetscFunctionReturn(0);
502d763cef2SBarry Smith }
503d763cef2SBarry Smith 
5044a2ae208SSatish Balay #undef __FUNCT__
5054a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
506d763cef2SBarry Smith /*@C
507d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
508d763cef2SBarry Smith     boundary conditions for the function F.
509d763cef2SBarry Smith 
510d763cef2SBarry Smith     Collective on TS
511d763cef2SBarry Smith 
512d763cef2SBarry Smith     Input Parameters:
513d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
514d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
515d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
516d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
517d763cef2SBarry Smith 
518d763cef2SBarry Smith     Calling sequence of func:
51987828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
520d763cef2SBarry Smith 
521d763cef2SBarry Smith +   t - current timestep
522d763cef2SBarry Smith .   F - function vector
523d763cef2SBarry Smith -   ctx - [optional] user-defined function context
524d763cef2SBarry Smith 
525d763cef2SBarry Smith     Level: intermediate
526d763cef2SBarry Smith 
527d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
528d763cef2SBarry Smith @*/
52987828ca2SBarry Smith int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,PetscReal,Vec,void*),void *ctx)
530d763cef2SBarry Smith {
531d763cef2SBarry Smith   PetscFunctionBegin;
532d763cef2SBarry Smith 
533d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
534d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
53529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
536d763cef2SBarry Smith   }
537000e7ae3SMatthew Knepley   ts->ops->rhsbc = f;
538d763cef2SBarry Smith   ts->bcP        = ctx;
539d763cef2SBarry Smith   PetscFunctionReturn(0);
540d763cef2SBarry Smith }
541d763cef2SBarry Smith 
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "TSView"
5447e2c5f70SBarry Smith /*@C
545d763cef2SBarry Smith     TSView - Prints the TS data structure.
546d763cef2SBarry Smith 
5474c49b128SBarry Smith     Collective on TS
548d763cef2SBarry Smith 
549d763cef2SBarry Smith     Input Parameters:
550d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
551d763cef2SBarry Smith -   viewer - visualization context
552d763cef2SBarry Smith 
553d763cef2SBarry Smith     Options Database Key:
554d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
555d763cef2SBarry Smith 
556d763cef2SBarry Smith     Notes:
557d763cef2SBarry Smith     The available visualization contexts include
558b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
559b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
560d763cef2SBarry Smith          output where only the first processor opens
561d763cef2SBarry Smith          the file.  All other processors send their
562d763cef2SBarry Smith          data to the first processor to print.
563d763cef2SBarry Smith 
564d763cef2SBarry Smith     The user can open an alternative visualization context with
565b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
566d763cef2SBarry Smith 
567d763cef2SBarry Smith     Level: beginner
568d763cef2SBarry Smith 
569d763cef2SBarry Smith .keywords: TS, timestep, view
570d763cef2SBarry Smith 
571b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
572d763cef2SBarry Smith @*/
573b0a32e0cSBarry Smith int TSView(TS ts,PetscViewer viewer)
574d763cef2SBarry Smith {
575d763cef2SBarry Smith   int        ierr;
576454a90a3SBarry Smith   char       *type;
5776831982aSBarry Smith   PetscTruth isascii,isstring;
578d763cef2SBarry Smith 
579d763cef2SBarry Smith   PetscFunctionBegin;
580d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
581b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
582b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
5836831982aSBarry Smith   PetscCheckSameComm(ts,viewer);
584fd16b177SBarry Smith 
585b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
586b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
5870f5bd95cSBarry Smith   if (isascii) {
588b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
589454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
590454a90a3SBarry Smith     if (type) {
591b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
592184914b5SBarry Smith     } else {
593b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
594184914b5SBarry Smith     }
595000e7ae3SMatthew Knepley     if (ts->ops->view) {
596b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
597000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
598b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
599d763cef2SBarry Smith     }
600b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
601b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
602d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
603b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
604d763cef2SBarry Smith     }
605b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
6060f5bd95cSBarry Smith   } else if (isstring) {
607454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
608b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
609d763cef2SBarry Smith   }
610b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
61194b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
612d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
613b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
614d763cef2SBarry Smith   PetscFunctionReturn(0);
615d763cef2SBarry Smith }
616d763cef2SBarry Smith 
617d763cef2SBarry Smith 
6184a2ae208SSatish Balay #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
620d763cef2SBarry Smith /*@C
621d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
622d763cef2SBarry Smith    the timesteppers.
623d763cef2SBarry Smith 
624d763cef2SBarry Smith    Collective on TS
625d763cef2SBarry Smith 
626d763cef2SBarry Smith    Input Parameters:
627d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
628d763cef2SBarry Smith -  usrP - optional user context
629d763cef2SBarry Smith 
630d763cef2SBarry Smith    Level: intermediate
631d763cef2SBarry Smith 
632d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
633d763cef2SBarry Smith 
634d763cef2SBarry Smith .seealso: TSGetApplicationContext()
635d763cef2SBarry Smith @*/
636d763cef2SBarry Smith int TSSetApplicationContext(TS ts,void *usrP)
637d763cef2SBarry Smith {
638d763cef2SBarry Smith   PetscFunctionBegin;
639d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
640d763cef2SBarry Smith   ts->user = usrP;
641d763cef2SBarry Smith   PetscFunctionReturn(0);
642d763cef2SBarry Smith }
643d763cef2SBarry Smith 
6444a2ae208SSatish Balay #undef __FUNCT__
6454a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
646d763cef2SBarry Smith /*@C
647d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
648d763cef2SBarry Smith     timestepper.
649d763cef2SBarry Smith 
650d763cef2SBarry Smith     Not Collective
651d763cef2SBarry Smith 
652d763cef2SBarry Smith     Input Parameter:
653d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
654d763cef2SBarry Smith 
655d763cef2SBarry Smith     Output Parameter:
656d763cef2SBarry Smith .   usrP - user context
657d763cef2SBarry Smith 
658d763cef2SBarry Smith     Level: intermediate
659d763cef2SBarry Smith 
660d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
661d763cef2SBarry Smith 
662d763cef2SBarry Smith .seealso: TSSetApplicationContext()
663d763cef2SBarry Smith @*/
664d763cef2SBarry Smith int TSGetApplicationContext(TS ts,void **usrP)
665d763cef2SBarry Smith {
666d763cef2SBarry Smith   PetscFunctionBegin;
667d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
668d763cef2SBarry Smith   *usrP = ts->user;
669d763cef2SBarry Smith   PetscFunctionReturn(0);
670d763cef2SBarry Smith }
671d763cef2SBarry Smith 
6724a2ae208SSatish Balay #undef __FUNCT__
6734a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
674d763cef2SBarry Smith /*@
675d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
676d763cef2SBarry Smith 
677d763cef2SBarry Smith    Not Collective
678d763cef2SBarry Smith 
679d763cef2SBarry Smith    Input Parameter:
680d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
681d763cef2SBarry Smith 
682d763cef2SBarry Smith    Output Parameter:
683d763cef2SBarry Smith .  iter - number steps so far
684d763cef2SBarry Smith 
685d763cef2SBarry Smith    Level: intermediate
686d763cef2SBarry Smith 
687d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
688d763cef2SBarry Smith @*/
689d763cef2SBarry Smith int TSGetTimeStepNumber(TS ts,int* iter)
690d763cef2SBarry Smith {
691d763cef2SBarry Smith   PetscFunctionBegin;
692d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
693d763cef2SBarry Smith   *iter = ts->steps;
694d763cef2SBarry Smith   PetscFunctionReturn(0);
695d763cef2SBarry Smith }
696d763cef2SBarry Smith 
6974a2ae208SSatish Balay #undef __FUNCT__
6984a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
699d763cef2SBarry Smith /*@
700d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
701d763cef2SBarry Smith    as well as the initial time.
702d763cef2SBarry Smith 
703d763cef2SBarry Smith    Collective on TS
704d763cef2SBarry Smith 
705d763cef2SBarry Smith    Input Parameters:
706d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
707d763cef2SBarry Smith .  initial_time - the initial time
708d763cef2SBarry Smith -  time_step - the size of the timestep
709d763cef2SBarry Smith 
710d763cef2SBarry Smith    Level: intermediate
711d763cef2SBarry Smith 
712d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
713d763cef2SBarry Smith 
714d763cef2SBarry Smith .keywords: TS, set, initial, timestep
715d763cef2SBarry Smith @*/
71687828ca2SBarry Smith int TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
717d763cef2SBarry Smith {
718d763cef2SBarry Smith   PetscFunctionBegin;
719d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
720d763cef2SBarry Smith   ts->time_step         = time_step;
721d763cef2SBarry Smith   ts->initial_time_step = time_step;
722d763cef2SBarry Smith   ts->ptime             = initial_time;
723d763cef2SBarry Smith   PetscFunctionReturn(0);
724d763cef2SBarry Smith }
725d763cef2SBarry Smith 
7264a2ae208SSatish Balay #undef __FUNCT__
7274a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
728d763cef2SBarry Smith /*@
729d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
730d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
731d763cef2SBarry Smith 
732d763cef2SBarry Smith    Collective on TS
733d763cef2SBarry Smith 
734d763cef2SBarry Smith    Input Parameters:
735d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
736d763cef2SBarry Smith -  time_step - the size of the timestep
737d763cef2SBarry Smith 
738d763cef2SBarry Smith    Level: intermediate
739d763cef2SBarry Smith 
740d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
741d763cef2SBarry Smith 
742d763cef2SBarry Smith .keywords: TS, set, timestep
743d763cef2SBarry Smith @*/
74487828ca2SBarry Smith int TSSetTimeStep(TS ts,PetscReal time_step)
745d763cef2SBarry Smith {
746d763cef2SBarry Smith   PetscFunctionBegin;
747d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
748d763cef2SBarry Smith   ts->time_step = time_step;
749d763cef2SBarry Smith   PetscFunctionReturn(0);
750d763cef2SBarry Smith }
751d763cef2SBarry Smith 
7524a2ae208SSatish Balay #undef __FUNCT__
7534a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
754d763cef2SBarry Smith /*@
755d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
756d763cef2SBarry Smith 
757d763cef2SBarry Smith    Not Collective
758d763cef2SBarry Smith 
759d763cef2SBarry Smith    Input Parameter:
760d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
761d763cef2SBarry Smith 
762d763cef2SBarry Smith    Output Parameter:
763d763cef2SBarry Smith .  dt - the current timestep size
764d763cef2SBarry Smith 
765d763cef2SBarry Smith    Level: intermediate
766d763cef2SBarry Smith 
767d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
768d763cef2SBarry Smith 
769d763cef2SBarry Smith .keywords: TS, get, timestep
770d763cef2SBarry Smith @*/
77187828ca2SBarry Smith int TSGetTimeStep(TS ts,PetscReal* dt)
772d763cef2SBarry Smith {
773d763cef2SBarry Smith   PetscFunctionBegin;
774d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
775d763cef2SBarry Smith   *dt = ts->time_step;
776d763cef2SBarry Smith   PetscFunctionReturn(0);
777d763cef2SBarry Smith }
778d763cef2SBarry Smith 
7794a2ae208SSatish Balay #undef __FUNCT__
7804a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
781d763cef2SBarry Smith /*@C
782d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
783d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
784d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
785d763cef2SBarry Smith    the solution at the next timestep has been calculated.
786d763cef2SBarry Smith 
787d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
788d763cef2SBarry Smith 
789d763cef2SBarry Smith    Input Parameter:
790d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
791d763cef2SBarry Smith 
792d763cef2SBarry Smith    Output Parameter:
793d763cef2SBarry Smith .  v - the vector containing the solution
794d763cef2SBarry Smith 
795d763cef2SBarry Smith    Level: intermediate
796d763cef2SBarry Smith 
797d763cef2SBarry Smith .seealso: TSGetTimeStep()
798d763cef2SBarry Smith 
799d763cef2SBarry Smith .keywords: TS, timestep, get, solution
800d763cef2SBarry Smith @*/
801d763cef2SBarry Smith int TSGetSolution(TS ts,Vec *v)
802d763cef2SBarry Smith {
803d763cef2SBarry Smith   PetscFunctionBegin;
804d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
805d763cef2SBarry Smith   *v = ts->vec_sol_always;
806d763cef2SBarry Smith   PetscFunctionReturn(0);
807d763cef2SBarry Smith }
808d763cef2SBarry Smith 
809bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
8104a2ae208SSatish Balay #undef __FUNCT__
811bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
812d763cef2SBarry Smith /*@C
813bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
814d763cef2SBarry Smith 
815bdad233fSMatthew Knepley   Not collective
816d763cef2SBarry Smith 
817bdad233fSMatthew Knepley   Input Parameters:
818bdad233fSMatthew Knepley + ts   - The TS
819bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
820d763cef2SBarry Smith .vb
821d763cef2SBarry Smith          U_t = A U
822d763cef2SBarry Smith          U_t = A(t) U
823d763cef2SBarry Smith          U_t = F(t,U)
824d763cef2SBarry Smith .ve
825d763cef2SBarry Smith 
826d763cef2SBarry Smith    Level: beginner
827d763cef2SBarry Smith 
828bdad233fSMatthew Knepley .keywords: TS, problem type
829bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
830d763cef2SBarry Smith @*/
831bdad233fSMatthew Knepley int TSSetProblemType(TS ts, TSProblemType type) {
832d763cef2SBarry Smith   PetscFunctionBegin;
833bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
834bdad233fSMatthew Knepley   ts->problem_type = type;
835d763cef2SBarry Smith   PetscFunctionReturn(0);
836d763cef2SBarry Smith }
837d763cef2SBarry Smith 
838bdad233fSMatthew Knepley #undef __FUNCT__
839bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
840bdad233fSMatthew Knepley /*@C
841bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
842bdad233fSMatthew Knepley 
843bdad233fSMatthew Knepley   Not collective
844bdad233fSMatthew Knepley 
845bdad233fSMatthew Knepley   Input Parameter:
846bdad233fSMatthew Knepley . ts   - The TS
847bdad233fSMatthew Knepley 
848bdad233fSMatthew Knepley   Output Parameter:
849bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
850bdad233fSMatthew Knepley .vb
851bdad233fSMatthew Knepley          U_t = A U
852bdad233fSMatthew Knepley          U_t = A(t) U
853bdad233fSMatthew Knepley          U_t = F(t,U)
854bdad233fSMatthew Knepley .ve
855bdad233fSMatthew Knepley 
856bdad233fSMatthew Knepley    Level: beginner
857bdad233fSMatthew Knepley 
858bdad233fSMatthew Knepley .keywords: TS, problem type
859bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
860bdad233fSMatthew Knepley @*/
861bdad233fSMatthew Knepley int TSGetProblemType(TS ts, TSProblemType *type) {
862bdad233fSMatthew Knepley   PetscFunctionBegin;
863bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
864bdad233fSMatthew Knepley   PetscValidPointer(type);
865bdad233fSMatthew Knepley   *type = ts->problem_type;
866bdad233fSMatthew Knepley   PetscFunctionReturn(0);
867bdad233fSMatthew Knepley }
868d763cef2SBarry Smith 
8694a2ae208SSatish Balay #undef __FUNCT__
8704a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
871d763cef2SBarry Smith /*@
872d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
873d763cef2SBarry Smith    of a timestepper.
874d763cef2SBarry Smith 
875d763cef2SBarry Smith    Collective on TS
876d763cef2SBarry Smith 
877d763cef2SBarry Smith    Input Parameter:
878d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
879d763cef2SBarry Smith 
880d763cef2SBarry Smith    Notes:
881d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
882d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
883d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
884d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
885d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
886d763cef2SBarry Smith 
887d763cef2SBarry Smith    Level: advanced
888d763cef2SBarry Smith 
889d763cef2SBarry Smith .keywords: TS, timestep, setup
890d763cef2SBarry Smith 
891d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
892d763cef2SBarry Smith @*/
893d763cef2SBarry Smith int TSSetUp(TS ts)
894d763cef2SBarry Smith {
895d763cef2SBarry Smith   int ierr;
896d763cef2SBarry Smith 
897d763cef2SBarry Smith   PetscFunctionBegin;
898d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
89929bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
900d763cef2SBarry Smith   if (!ts->type_name) {
901d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
902d763cef2SBarry Smith   }
903000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
904d763cef2SBarry Smith   ts->setupcalled = 1;
905d763cef2SBarry Smith   PetscFunctionReturn(0);
906d763cef2SBarry Smith }
907d763cef2SBarry Smith 
9084a2ae208SSatish Balay #undef __FUNCT__
9094a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
910d763cef2SBarry Smith /*@C
911d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
912d763cef2SBarry Smith    with TSCreate().
913d763cef2SBarry Smith 
914d763cef2SBarry Smith    Collective on TS
915d763cef2SBarry Smith 
916d763cef2SBarry Smith    Input Parameter:
917d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
918d763cef2SBarry Smith 
919d763cef2SBarry Smith    Level: beginner
920d763cef2SBarry Smith 
921d763cef2SBarry Smith .keywords: TS, timestepper, destroy
922d763cef2SBarry Smith 
923d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
924d763cef2SBarry Smith @*/
925d763cef2SBarry Smith int TSDestroy(TS ts)
926d763cef2SBarry Smith {
927329f5518SBarry Smith   int ierr,i;
928d763cef2SBarry Smith 
929d763cef2SBarry Smith   PetscFunctionBegin;
930d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
931d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
932d763cef2SBarry Smith 
933be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
9340f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
935be0abb6dSBarry Smith 
93694b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
937d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
938000e7ae3SMatthew Knepley   ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);
939329f5518SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
940329f5518SBarry Smith     if (ts->mdestroy[i]) {
941329f5518SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
942329f5518SBarry Smith     }
943329f5518SBarry Smith   }
944b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)ts);
945d763cef2SBarry Smith   PetscHeaderDestroy((PetscObject)ts);
946d763cef2SBarry Smith   PetscFunctionReturn(0);
947d763cef2SBarry Smith }
948d763cef2SBarry Smith 
9494a2ae208SSatish Balay #undef __FUNCT__
9504a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
951d763cef2SBarry Smith /*@C
952d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
953d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
954d763cef2SBarry Smith 
955d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
956d763cef2SBarry Smith 
957d763cef2SBarry Smith    Input Parameter:
958d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
959d763cef2SBarry Smith 
960d763cef2SBarry Smith    Output Parameter:
961d763cef2SBarry Smith .  snes - the nonlinear solver context
962d763cef2SBarry Smith 
963d763cef2SBarry Smith    Notes:
964d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
965d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
96694b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
967d763cef2SBarry Smith 
968d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
969d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
970d763cef2SBarry Smith 
971d763cef2SBarry Smith    Level: beginner
972d763cef2SBarry Smith 
973d763cef2SBarry Smith .keywords: timestep, get, SNES
974d763cef2SBarry Smith @*/
975d763cef2SBarry Smith int TSGetSNES(TS ts,SNES *snes)
976d763cef2SBarry Smith {
977d763cef2SBarry Smith   PetscFunctionBegin;
978d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
97994b7f48cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
980d763cef2SBarry Smith   *snes = ts->snes;
981d763cef2SBarry Smith   PetscFunctionReturn(0);
982d763cef2SBarry Smith }
983d763cef2SBarry Smith 
9844a2ae208SSatish Balay #undef __FUNCT__
98594b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
986d763cef2SBarry Smith /*@C
98794b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
988d763cef2SBarry Smith    a TS (timestepper) context.
989d763cef2SBarry Smith 
99094b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
991d763cef2SBarry Smith 
992d763cef2SBarry Smith    Input Parameter:
993d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
994d763cef2SBarry Smith 
995d763cef2SBarry Smith    Output Parameter:
99694b7f48cSBarry Smith .  ksp - the nonlinear solver context
997d763cef2SBarry Smith 
998d763cef2SBarry Smith    Notes:
99994b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1000d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1001d763cef2SBarry Smith    KSP and PC contexts as well.
1002d763cef2SBarry Smith 
100394b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
100494b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1005d763cef2SBarry Smith 
1006d763cef2SBarry Smith    Level: beginner
1007d763cef2SBarry Smith 
100894b7f48cSBarry Smith .keywords: timestep, get, KSP
1009d763cef2SBarry Smith @*/
101094b7f48cSBarry Smith int TSGetKSP(TS ts,KSP *ksp)
1011d763cef2SBarry Smith {
1012d763cef2SBarry Smith   PetscFunctionBegin;
1013d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
101429bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
101594b7f48cSBarry Smith   *ksp = ts->ksp;
1016d763cef2SBarry Smith   PetscFunctionReturn(0);
1017d763cef2SBarry Smith }
1018d763cef2SBarry Smith 
1019d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1020d763cef2SBarry Smith 
10214a2ae208SSatish Balay #undef __FUNCT__
1022adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1023adb62b0dSMatthew Knepley /*@
1024adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1025adb62b0dSMatthew Knepley    maximum time for iteration.
1026adb62b0dSMatthew Knepley 
1027adb62b0dSMatthew Knepley    Collective on TS
1028adb62b0dSMatthew Knepley 
1029adb62b0dSMatthew Knepley    Input Parameters:
1030adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1031adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1032adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1033adb62b0dSMatthew Knepley 
1034adb62b0dSMatthew Knepley    Level: intermediate
1035adb62b0dSMatthew Knepley 
1036adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1037adb62b0dSMatthew Knepley @*/
1038adb62b0dSMatthew Knepley int TSGetDuration(TS ts, int *maxsteps, PetscReal *maxtime)
1039adb62b0dSMatthew Knepley {
1040adb62b0dSMatthew Knepley   PetscFunctionBegin;
1041adb62b0dSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1042adb62b0dSMatthew Knepley   if (maxsteps != PETSC_NULL) {
1043adb62b0dSMatthew Knepley     PetscValidIntPointer(maxsteps);
1044adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1045adb62b0dSMatthew Knepley   }
1046adb62b0dSMatthew Knepley   if (maxtime  != PETSC_NULL) {
1047adb62b0dSMatthew Knepley     PetscValidScalarPointer(maxtime);
1048adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1049adb62b0dSMatthew Knepley   }
1050adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1051adb62b0dSMatthew Knepley }
1052adb62b0dSMatthew Knepley 
1053adb62b0dSMatthew Knepley #undef __FUNCT__
10544a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1055d763cef2SBarry Smith /*@
1056d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1057d763cef2SBarry Smith    maximum time for iteration.
1058d763cef2SBarry Smith 
1059d763cef2SBarry Smith    Collective on TS
1060d763cef2SBarry Smith 
1061d763cef2SBarry Smith    Input Parameters:
1062d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1063d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1064d763cef2SBarry Smith -  maxtime - final time to iterate to
1065d763cef2SBarry Smith 
1066d763cef2SBarry Smith    Options Database Keys:
1067d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1068d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1069d763cef2SBarry Smith 
1070d763cef2SBarry Smith    Notes:
1071d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1072d763cef2SBarry Smith 
1073d763cef2SBarry Smith    Level: intermediate
1074d763cef2SBarry Smith 
1075d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1076d763cef2SBarry Smith @*/
107787828ca2SBarry Smith int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime)
1078d763cef2SBarry Smith {
1079d763cef2SBarry Smith   PetscFunctionBegin;
1080d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1081d763cef2SBarry Smith   ts->max_steps = maxsteps;
1082d763cef2SBarry Smith   ts->max_time  = maxtime;
1083d763cef2SBarry Smith   PetscFunctionReturn(0);
1084d763cef2SBarry Smith }
1085d763cef2SBarry Smith 
10864a2ae208SSatish Balay #undef __FUNCT__
10874a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1088d763cef2SBarry Smith /*@
1089d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1090d763cef2SBarry Smith    for use by the TS routines.
1091d763cef2SBarry Smith 
1092d763cef2SBarry Smith    Collective on TS and Vec
1093d763cef2SBarry Smith 
1094d763cef2SBarry Smith    Input Parameters:
1095d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1096d763cef2SBarry Smith -  x - the solution vector
1097d763cef2SBarry Smith 
1098d763cef2SBarry Smith    Level: beginner
1099d763cef2SBarry Smith 
1100d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1101d763cef2SBarry Smith @*/
1102d763cef2SBarry Smith int TSSetSolution(TS ts,Vec x)
1103d763cef2SBarry Smith {
1104d763cef2SBarry Smith   PetscFunctionBegin;
1105d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1106d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1107d763cef2SBarry Smith   PetscFunctionReturn(0);
1108d763cef2SBarry Smith }
1109d763cef2SBarry Smith 
1110e74ef692SMatthew Knepley #undef __FUNCT__
1111e74ef692SMatthew Knepley #define __FUNCT__ "TSSetRhsBC"
1112000e7ae3SMatthew Knepley /*@
1113000e7ae3SMatthew Knepley   TSSetRhsBC - Sets the function which applies boundary conditions
1114000e7ae3SMatthew Knepley   to the Rhs of each system.
1115000e7ae3SMatthew Knepley 
1116000e7ae3SMatthew Knepley   Collective on TS
1117000e7ae3SMatthew Knepley 
1118000e7ae3SMatthew Knepley   Input Parameters:
1119000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1120000e7ae3SMatthew Knepley - func - The function
1121000e7ae3SMatthew Knepley 
1122000e7ae3SMatthew Knepley   Calling sequence of func:
1123000e7ae3SMatthew Knepley . func (TS ts, Vec rhs, void *ctx);
1124000e7ae3SMatthew Knepley 
1125000e7ae3SMatthew Knepley + rhs - The current rhs vector
1126000e7ae3SMatthew Knepley - ctx - The user-context
1127000e7ae3SMatthew Knepley 
1128000e7ae3SMatthew Knepley   Level: intermediate
1129000e7ae3SMatthew Knepley 
1130000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1131000e7ae3SMatthew Knepley @*/
1132000e7ae3SMatthew Knepley int TSSetRhsBC(TS ts, int (*func)(TS, Vec, void *))
1133000e7ae3SMatthew Knepley {
1134000e7ae3SMatthew Knepley   PetscFunctionBegin;
1135000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1136000e7ae3SMatthew Knepley   ts->ops->applyrhsbc = func;
1137000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1138000e7ae3SMatthew Knepley }
1139000e7ae3SMatthew Knepley 
1140e74ef692SMatthew Knepley #undef __FUNCT__
1141e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultRhsBC"
1142000e7ae3SMatthew Knepley /*@
1143000e7ae3SMatthew Knepley   TSDefaultRhsBC - The default boundary condition function which does nothing.
1144000e7ae3SMatthew Knepley 
1145000e7ae3SMatthew Knepley   Collective on TS
1146000e7ae3SMatthew Knepley 
1147000e7ae3SMatthew Knepley   Input Parameters:
1148000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1149000e7ae3SMatthew Knepley . rhs - The Rhs
1150000e7ae3SMatthew Knepley - ctx - The user-context
1151000e7ae3SMatthew Knepley 
1152000e7ae3SMatthew Knepley   Level: developer
1153000e7ae3SMatthew Knepley 
1154000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1155000e7ae3SMatthew Knepley @*/
1156000e7ae3SMatthew Knepley int TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1157000e7ae3SMatthew Knepley {
1158000e7ae3SMatthew Knepley   PetscFunctionBegin;
1159000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1160000e7ae3SMatthew Knepley }
1161000e7ae3SMatthew Knepley 
1162e74ef692SMatthew Knepley #undef __FUNCT__
1163e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSystemMatrixBC"
1164000e7ae3SMatthew Knepley /*@
1165000e7ae3SMatthew Knepley   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1166000e7ae3SMatthew Knepley   to the system matrix and preconditioner of each system.
1167000e7ae3SMatthew Knepley 
1168000e7ae3SMatthew Knepley   Collective on TS
1169000e7ae3SMatthew Knepley 
1170000e7ae3SMatthew Knepley   Input Parameters:
1171000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1172000e7ae3SMatthew Knepley - func - The function
1173000e7ae3SMatthew Knepley 
1174000e7ae3SMatthew Knepley   Calling sequence of func:
1175000e7ae3SMatthew Knepley . func (TS ts, Mat A, Mat B, void *ctx);
1176000e7ae3SMatthew Knepley 
1177000e7ae3SMatthew Knepley + A   - The current system matrix
1178000e7ae3SMatthew Knepley . B   - The current preconditioner
1179000e7ae3SMatthew Knepley - ctx - The user-context
1180000e7ae3SMatthew Knepley 
1181000e7ae3SMatthew Knepley   Level: intermediate
1182000e7ae3SMatthew Knepley 
1183000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1184000e7ae3SMatthew Knepley @*/
1185000e7ae3SMatthew Knepley int TSSetSystemMatrixBC(TS ts, int (*func)(TS, Mat, Mat, void *))
1186000e7ae3SMatthew Knepley {
1187000e7ae3SMatthew Knepley   PetscFunctionBegin;
1188000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1189000e7ae3SMatthew Knepley   ts->ops->applymatrixbc = func;
1190000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1191000e7ae3SMatthew Knepley }
1192000e7ae3SMatthew Knepley 
1193e74ef692SMatthew Knepley #undef __FUNCT__
1194e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSystemMatrixBC"
1195000e7ae3SMatthew Knepley /*@
1196000e7ae3SMatthew Knepley   TSDefaultSystemMatrixBC - The default boundary condition function which
1197000e7ae3SMatthew Knepley   does nothing.
1198000e7ae3SMatthew Knepley 
1199000e7ae3SMatthew Knepley   Collective on TS
1200000e7ae3SMatthew Knepley 
1201000e7ae3SMatthew Knepley   Input Parameters:
1202000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1203000e7ae3SMatthew Knepley . A   - The system matrix
1204000e7ae3SMatthew Knepley . B   - The preconditioner
1205000e7ae3SMatthew Knepley - ctx - The user-context
1206000e7ae3SMatthew Knepley 
1207000e7ae3SMatthew Knepley   Level: developer
1208000e7ae3SMatthew Knepley 
1209000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1210000e7ae3SMatthew Knepley @*/
1211000e7ae3SMatthew Knepley int TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1212000e7ae3SMatthew Knepley {
1213000e7ae3SMatthew Knepley   PetscFunctionBegin;
1214000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1215000e7ae3SMatthew Knepley }
1216000e7ae3SMatthew Knepley 
1217e74ef692SMatthew Knepley #undef __FUNCT__
1218e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSolutionBC"
1219000e7ae3SMatthew Knepley /*@
1220000e7ae3SMatthew Knepley   TSSetSolutionBC - Sets the function which applies boundary conditions
1221000e7ae3SMatthew Knepley   to the solution of each system. This is necessary in nonlinear systems
1222000e7ae3SMatthew Knepley   which time dependent boundary conditions.
1223000e7ae3SMatthew Knepley 
1224000e7ae3SMatthew Knepley   Collective on TS
1225000e7ae3SMatthew Knepley 
1226000e7ae3SMatthew Knepley   Input Parameters:
1227000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1228000e7ae3SMatthew Knepley - func - The function
1229000e7ae3SMatthew Knepley 
1230000e7ae3SMatthew Knepley   Calling sequence of func:
1231000e7ae3SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
1232000e7ae3SMatthew Knepley 
1233000e7ae3SMatthew Knepley + sol - The current solution vector
1234000e7ae3SMatthew Knepley - ctx - The user-context
1235000e7ae3SMatthew Knepley 
1236000e7ae3SMatthew Knepley   Level: intermediate
1237000e7ae3SMatthew Knepley 
1238000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1239000e7ae3SMatthew Knepley @*/
1240000e7ae3SMatthew Knepley int TSSetSolutionBC(TS ts, int (*func)(TS, Vec, void *))
1241000e7ae3SMatthew Knepley {
1242000e7ae3SMatthew Knepley   PetscFunctionBegin;
1243000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1244000e7ae3SMatthew Knepley   ts->ops->applysolbc = func;
1245000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1246000e7ae3SMatthew Knepley }
1247000e7ae3SMatthew Knepley 
1248e74ef692SMatthew Knepley #undef __FUNCT__
1249e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSolutionBC"
1250000e7ae3SMatthew Knepley /*@
1251000e7ae3SMatthew Knepley   TSDefaultSolutionBC - The default boundary condition function which
1252000e7ae3SMatthew Knepley   does nothing.
1253000e7ae3SMatthew Knepley 
1254000e7ae3SMatthew Knepley   Collective on TS
1255000e7ae3SMatthew Knepley 
1256000e7ae3SMatthew Knepley   Input Parameters:
1257000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1258000e7ae3SMatthew Knepley . sol - The solution
1259000e7ae3SMatthew Knepley - ctx - The user-context
1260000e7ae3SMatthew Knepley 
1261000e7ae3SMatthew Knepley   Level: developer
1262000e7ae3SMatthew Knepley 
1263000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1264000e7ae3SMatthew Knepley @*/
1265000e7ae3SMatthew Knepley int TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1266000e7ae3SMatthew Knepley {
1267000e7ae3SMatthew Knepley   PetscFunctionBegin;
1268000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1269000e7ae3SMatthew Knepley }
1270000e7ae3SMatthew Knepley 
1271e74ef692SMatthew Knepley #undef __FUNCT__
1272e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1273000e7ae3SMatthew Knepley /*@
1274000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1275000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1276000e7ae3SMatthew Knepley 
1277000e7ae3SMatthew Knepley   Collective on TS
1278000e7ae3SMatthew Knepley 
1279000e7ae3SMatthew Knepley   Input Parameters:
1280000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1281000e7ae3SMatthew Knepley - func - The function
1282000e7ae3SMatthew Knepley 
1283000e7ae3SMatthew Knepley   Calling sequence of func:
1284000e7ae3SMatthew Knepley . func (TS ts);
1285000e7ae3SMatthew Knepley 
1286000e7ae3SMatthew Knepley   Level: intermediate
1287000e7ae3SMatthew Knepley 
1288000e7ae3SMatthew Knepley .keywords: TS, timestep
1289000e7ae3SMatthew Knepley @*/
1290000e7ae3SMatthew Knepley int TSSetPreStep(TS ts, int (*func)(TS))
1291000e7ae3SMatthew Knepley {
1292000e7ae3SMatthew Knepley   PetscFunctionBegin;
1293000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1294000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1295000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1296000e7ae3SMatthew Knepley }
1297000e7ae3SMatthew Knepley 
1298e74ef692SMatthew Knepley #undef __FUNCT__
1299e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1300000e7ae3SMatthew Knepley /*@
1301000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1302000e7ae3SMatthew Knepley 
1303000e7ae3SMatthew Knepley   Collective on TS
1304000e7ae3SMatthew Knepley 
1305000e7ae3SMatthew Knepley   Input Parameters:
1306000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1307000e7ae3SMatthew Knepley 
1308000e7ae3SMatthew Knepley   Level: developer
1309000e7ae3SMatthew Knepley 
1310000e7ae3SMatthew Knepley .keywords: TS, timestep
1311000e7ae3SMatthew Knepley @*/
1312000e7ae3SMatthew Knepley int TSDefaultPreStep(TS ts)
1313000e7ae3SMatthew Knepley {
1314000e7ae3SMatthew Knepley   PetscFunctionBegin;
1315000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1316000e7ae3SMatthew Knepley }
1317000e7ae3SMatthew Knepley 
1318e74ef692SMatthew Knepley #undef __FUNCT__
1319e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1320000e7ae3SMatthew Knepley /*@
1321000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1322000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1323000e7ae3SMatthew Knepley   the time step.
1324000e7ae3SMatthew Knepley 
1325000e7ae3SMatthew Knepley   Collective on TS
1326000e7ae3SMatthew Knepley 
1327000e7ae3SMatthew Knepley   Input Parameters:
1328000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1329000e7ae3SMatthew Knepley - func - The function
1330000e7ae3SMatthew Knepley 
1331000e7ae3SMatthew Knepley   Calling sequence of func:
1332000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1333000e7ae3SMatthew Knepley 
1334000e7ae3SMatthew Knepley + t   - The current time
1335000e7ae3SMatthew Knepley - dt  - The current time step
1336000e7ae3SMatthew Knepley 
1337000e7ae3SMatthew Knepley   Level: intermediate
1338000e7ae3SMatthew Knepley 
1339000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1340000e7ae3SMatthew Knepley @*/
134189cef881SMatthew Knepley int TSSetUpdate(TS ts, int (*func)(TS, PetscReal, PetscReal *))
1342000e7ae3SMatthew Knepley {
1343000e7ae3SMatthew Knepley   PetscFunctionBegin;
1344000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1345000e7ae3SMatthew Knepley   ts->ops->update = func;
1346000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1347000e7ae3SMatthew Knepley }
1348000e7ae3SMatthew Knepley 
1349e74ef692SMatthew Knepley #undef __FUNCT__
1350e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1351000e7ae3SMatthew Knepley /*@
1352000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1353000e7ae3SMatthew Knepley 
1354000e7ae3SMatthew Knepley   Collective on TS
1355000e7ae3SMatthew Knepley 
1356000e7ae3SMatthew Knepley   Input Parameters:
1357000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1358000e7ae3SMatthew Knepley - t   - The current time
1359000e7ae3SMatthew Knepley 
1360000e7ae3SMatthew Knepley   Output Parameters:
1361000e7ae3SMatthew Knepley . dt  - The current time step
1362000e7ae3SMatthew Knepley 
1363000e7ae3SMatthew Knepley   Level: developer
1364000e7ae3SMatthew Knepley 
1365000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1366000e7ae3SMatthew Knepley @*/
136789cef881SMatthew Knepley int TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1368000e7ae3SMatthew Knepley {
1369000e7ae3SMatthew Knepley   PetscFunctionBegin;
1370000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1371000e7ae3SMatthew Knepley }
1372000e7ae3SMatthew Knepley 
1373e74ef692SMatthew Knepley #undef __FUNCT__
1374e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1375000e7ae3SMatthew Knepley /*@
1376000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1377000e7ae3SMatthew Knepley   called once at the end of time stepping.
1378000e7ae3SMatthew Knepley 
1379000e7ae3SMatthew Knepley   Collective on TS
1380000e7ae3SMatthew Knepley 
1381000e7ae3SMatthew Knepley   Input Parameters:
1382000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1383000e7ae3SMatthew Knepley - func - The function
1384000e7ae3SMatthew Knepley 
1385000e7ae3SMatthew Knepley   Calling sequence of func:
1386000e7ae3SMatthew Knepley . func (TS ts);
1387000e7ae3SMatthew Knepley 
1388000e7ae3SMatthew Knepley   Level: intermediate
1389000e7ae3SMatthew Knepley 
1390000e7ae3SMatthew Knepley .keywords: TS, timestep
1391000e7ae3SMatthew Knepley @*/
1392000e7ae3SMatthew Knepley int TSSetPostStep(TS ts, int (*func)(TS))
1393000e7ae3SMatthew Knepley {
1394000e7ae3SMatthew Knepley   PetscFunctionBegin;
1395000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1396000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1397000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1398000e7ae3SMatthew Knepley }
1399000e7ae3SMatthew Knepley 
1400e74ef692SMatthew Knepley #undef __FUNCT__
1401e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1402000e7ae3SMatthew Knepley /*@
1403000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1404000e7ae3SMatthew Knepley 
1405000e7ae3SMatthew Knepley   Collective on TS
1406000e7ae3SMatthew Knepley 
1407000e7ae3SMatthew Knepley   Input Parameters:
1408000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1409000e7ae3SMatthew Knepley 
1410000e7ae3SMatthew Knepley   Level: developer
1411000e7ae3SMatthew Knepley 
1412000e7ae3SMatthew Knepley .keywords: TS, timestep
1413000e7ae3SMatthew Knepley @*/
1414000e7ae3SMatthew Knepley int TSDefaultPostStep(TS ts)
1415000e7ae3SMatthew Knepley {
1416000e7ae3SMatthew Knepley   PetscFunctionBegin;
1417000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1418000e7ae3SMatthew Knepley }
1419000e7ae3SMatthew Knepley 
1420d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1421d763cef2SBarry Smith 
14224a2ae208SSatish Balay #undef __FUNCT__
14234a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
1424d763cef2SBarry Smith /*@C
1425d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1426d763cef2SBarry Smith    timestep to display the iteration's  progress.
1427d763cef2SBarry Smith 
1428d763cef2SBarry Smith    Collective on TS
1429d763cef2SBarry Smith 
1430d763cef2SBarry Smith    Input Parameters:
1431d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1432d763cef2SBarry Smith .  func - monitoring routine
1433329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1434b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1435b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1436b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1437d763cef2SBarry Smith 
1438d763cef2SBarry Smith    Calling sequence of func:
143987828ca2SBarry Smith $    int func(TS ts,int steps,PetscReal time,Vec x,void *mctx)
1440d763cef2SBarry Smith 
1441d763cef2SBarry Smith +    ts - the TS context
1442d763cef2SBarry Smith .    steps - iteration number
14431f06c33eSBarry Smith .    time - current time
1444d763cef2SBarry Smith .    x - current iterate
1445d763cef2SBarry Smith -    mctx - [optional] monitoring context
1446d763cef2SBarry Smith 
1447d763cef2SBarry Smith    Notes:
1448d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1449d763cef2SBarry Smith    already has been loaded.
1450d763cef2SBarry Smith 
1451d763cef2SBarry Smith    Level: intermediate
1452d763cef2SBarry Smith 
1453d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1454d763cef2SBarry Smith 
1455d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
1456d763cef2SBarry Smith @*/
145787828ca2SBarry Smith int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*))
1458d763cef2SBarry Smith {
1459d763cef2SBarry Smith   PetscFunctionBegin;
1460d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1461d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
146229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1463d763cef2SBarry Smith   }
1464d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1465329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1466d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1467d763cef2SBarry Smith   PetscFunctionReturn(0);
1468d763cef2SBarry Smith }
1469d763cef2SBarry Smith 
14704a2ae208SSatish Balay #undef __FUNCT__
14714a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
1472d763cef2SBarry Smith /*@C
1473d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1474d763cef2SBarry Smith 
1475d763cef2SBarry Smith    Collective on TS
1476d763cef2SBarry Smith 
1477d763cef2SBarry Smith    Input Parameters:
1478d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1479d763cef2SBarry Smith 
1480d763cef2SBarry Smith    Notes:
1481d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1482d763cef2SBarry Smith 
1483d763cef2SBarry Smith    Level: intermediate
1484d763cef2SBarry Smith 
1485d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1486d763cef2SBarry Smith 
1487d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
1488d763cef2SBarry Smith @*/
1489d763cef2SBarry Smith int TSClearMonitor(TS ts)
1490d763cef2SBarry Smith {
1491d763cef2SBarry Smith   PetscFunctionBegin;
1492d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1493d763cef2SBarry Smith   ts->numbermonitors = 0;
1494d763cef2SBarry Smith   PetscFunctionReturn(0);
1495d763cef2SBarry Smith }
1496d763cef2SBarry Smith 
14974a2ae208SSatish Balay #undef __FUNCT__
14984a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
149987828ca2SBarry Smith int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
1500d763cef2SBarry Smith {
1501d132466eSBarry Smith   int ierr;
1502d132466eSBarry Smith 
1503d763cef2SBarry Smith   PetscFunctionBegin;
1504142b95e3SSatish Balay   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1505d763cef2SBarry Smith   PetscFunctionReturn(0);
1506d763cef2SBarry Smith }
1507d763cef2SBarry Smith 
15084a2ae208SSatish Balay #undef __FUNCT__
15094a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1510d763cef2SBarry Smith /*@
1511d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1512d763cef2SBarry Smith 
1513d763cef2SBarry Smith    Collective on TS
1514d763cef2SBarry Smith 
1515d763cef2SBarry Smith    Input Parameter:
1516d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1517d763cef2SBarry Smith 
1518d763cef2SBarry Smith    Output Parameters:
1519d763cef2SBarry Smith +  steps - number of iterations until termination
1520142b95e3SSatish Balay -  ptime - time until termination
1521d763cef2SBarry Smith 
1522d763cef2SBarry Smith    Level: beginner
1523d763cef2SBarry Smith 
1524d763cef2SBarry Smith .keywords: TS, timestep, solve
1525d763cef2SBarry Smith 
1526d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1527d763cef2SBarry Smith @*/
152887828ca2SBarry Smith int TSStep(TS ts,int *steps,PetscReal *ptime)
1529d763cef2SBarry Smith {
1530f1af5d2fSBarry Smith   int        ierr;
1531d763cef2SBarry Smith 
1532d763cef2SBarry Smith   PetscFunctionBegin;
1533d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1534d405a339SMatthew Knepley   if (!ts->setupcalled) {
1535d405a339SMatthew Knepley     ierr = TSSetUp(ts);                                                                                   CHKERRQ(ierr);
1536d405a339SMatthew Knepley   }
1537d405a339SMatthew Knepley 
1538d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);                                                        CHKERRQ(ierr);
1539d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);                                                                         CHKERRQ(ierr);
1540000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);                                                              CHKERRQ(ierr);
1541d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);                                                                        CHKERRQ(ierr);
1542d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);                                                          CHKERRQ(ierr);
1543d405a339SMatthew Knepley 
15444bb05414SBarry Smith   if (!PetscPreLoadingOn) {
15454bb05414SBarry Smith     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1546d405a339SMatthew Knepley   }
1547d763cef2SBarry Smith   PetscFunctionReturn(0);
1548d763cef2SBarry Smith }
1549d763cef2SBarry Smith 
15504a2ae208SSatish Balay #undef __FUNCT__
15514a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1552d763cef2SBarry Smith /*
1553d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1554d763cef2SBarry Smith */
155587828ca2SBarry Smith int TSMonitor(TS ts,int step,PetscReal ptime,Vec x)
1556d763cef2SBarry Smith {
1557d763cef2SBarry Smith   int i,ierr,n = ts->numbermonitors;
1558d763cef2SBarry Smith 
1559d763cef2SBarry Smith   PetscFunctionBegin;
1560d763cef2SBarry Smith   for (i=0; i<n; i++) {
1561142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1562d763cef2SBarry Smith   }
1563d763cef2SBarry Smith   PetscFunctionReturn(0);
1564d763cef2SBarry Smith }
1565d763cef2SBarry Smith 
1566d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1567d763cef2SBarry Smith 
15684a2ae208SSatish Balay #undef __FUNCT__
15694a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1570d763cef2SBarry Smith /*@C
1571d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1572d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1573d763cef2SBarry Smith 
1574d763cef2SBarry Smith    Collective on TS
1575d763cef2SBarry Smith 
1576d763cef2SBarry Smith    Input Parameters:
1577d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1578d763cef2SBarry Smith .  label - the title to put in the title bar
15797c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1580d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1581d763cef2SBarry Smith 
1582d763cef2SBarry Smith    Output Parameter:
1583d763cef2SBarry Smith .  draw - the drawing context
1584d763cef2SBarry Smith 
1585d763cef2SBarry Smith    Options Database Key:
1586d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1587d763cef2SBarry Smith 
1588d763cef2SBarry Smith    Notes:
1589b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1590d763cef2SBarry Smith 
1591d763cef2SBarry Smith    Level: intermediate
1592d763cef2SBarry Smith 
15937c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1594d763cef2SBarry Smith 
1595d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
15967c922b88SBarry Smith 
1597d763cef2SBarry Smith @*/
15981836bdbcSSatish Balay int TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1599d763cef2SBarry Smith {
1600b0a32e0cSBarry Smith   PetscDraw win;
1601d763cef2SBarry Smith   int       ierr;
1602d763cef2SBarry Smith 
1603d763cef2SBarry Smith   PetscFunctionBegin;
1604b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1605b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1606b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1607b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1608d763cef2SBarry Smith 
1609b0a32e0cSBarry Smith   PetscLogObjectParent(*draw,win);
1610d763cef2SBarry Smith   PetscFunctionReturn(0);
1611d763cef2SBarry Smith }
1612d763cef2SBarry Smith 
16134a2ae208SSatish Balay #undef __FUNCT__
16144a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
161587828ca2SBarry Smith int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx)
1616d763cef2SBarry Smith {
1617b0a32e0cSBarry Smith   PetscDrawLG lg = (PetscDrawLG) monctx;
161887828ca2SBarry Smith   PetscReal      x,y = ptime;
1619d763cef2SBarry Smith   int         ierr;
1620d763cef2SBarry Smith 
1621d763cef2SBarry Smith   PetscFunctionBegin;
16227c922b88SBarry Smith   if (!monctx) {
16237c922b88SBarry Smith     MPI_Comm    comm;
1624b0a32e0cSBarry Smith     PetscViewer viewer;
16257c922b88SBarry Smith 
16267c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1627b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1628b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
16297c922b88SBarry Smith   }
16307c922b88SBarry Smith 
1631b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
163287828ca2SBarry Smith   x = (PetscReal)n;
1633b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1634d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1635b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1636d763cef2SBarry Smith   }
1637d763cef2SBarry Smith   PetscFunctionReturn(0);
1638d763cef2SBarry Smith }
1639d763cef2SBarry Smith 
16404a2ae208SSatish Balay #undef __FUNCT__
16414a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1642d763cef2SBarry Smith /*@C
1643d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1644d763cef2SBarry Smith    with TSLGMonitorCreate().
1645d763cef2SBarry Smith 
1646b0a32e0cSBarry Smith    Collective on PetscDrawLG
1647d763cef2SBarry Smith 
1648d763cef2SBarry Smith    Input Parameter:
1649d763cef2SBarry Smith .  draw - the drawing context
1650d763cef2SBarry Smith 
1651d763cef2SBarry Smith    Level: intermediate
1652d763cef2SBarry Smith 
1653d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1654d763cef2SBarry Smith 
1655d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1656d763cef2SBarry Smith @*/
1657b0a32e0cSBarry Smith int TSLGMonitorDestroy(PetscDrawLG drawlg)
1658d763cef2SBarry Smith {
1659b0a32e0cSBarry Smith   PetscDraw draw;
1660d763cef2SBarry Smith   int       ierr;
1661d763cef2SBarry Smith 
1662d763cef2SBarry Smith   PetscFunctionBegin;
1663b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1664b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1665b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1666d763cef2SBarry Smith   PetscFunctionReturn(0);
1667d763cef2SBarry Smith }
1668d763cef2SBarry Smith 
16694a2ae208SSatish Balay #undef __FUNCT__
16704a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1671d763cef2SBarry Smith /*@
1672d763cef2SBarry Smith    TSGetTime - Gets the current time.
1673d763cef2SBarry Smith 
1674d763cef2SBarry Smith    Not Collective
1675d763cef2SBarry Smith 
1676d763cef2SBarry Smith    Input Parameter:
1677d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1678d763cef2SBarry Smith 
1679d763cef2SBarry Smith    Output Parameter:
1680d763cef2SBarry Smith .  t  - the current time
1681d763cef2SBarry Smith 
1682d763cef2SBarry Smith    Contributed by: Matthew Knepley
1683d763cef2SBarry Smith 
1684d763cef2SBarry Smith    Level: beginner
1685d763cef2SBarry Smith 
1686d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1687d763cef2SBarry Smith 
1688d763cef2SBarry Smith .keywords: TS, get, time
1689d763cef2SBarry Smith @*/
169087828ca2SBarry Smith int TSGetTime(TS ts,PetscReal* t)
1691d763cef2SBarry Smith {
1692d763cef2SBarry Smith   PetscFunctionBegin;
1693d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1694d763cef2SBarry Smith   *t = ts->ptime;
1695d763cef2SBarry Smith   PetscFunctionReturn(0);
1696d763cef2SBarry Smith }
1697d763cef2SBarry Smith 
16984a2ae208SSatish Balay #undef __FUNCT__
16994a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1700d763cef2SBarry Smith /*@C
1701d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1702d763cef2SBarry Smith    TS options in the database.
1703d763cef2SBarry Smith 
1704d763cef2SBarry Smith    Collective on TS
1705d763cef2SBarry Smith 
1706d763cef2SBarry Smith    Input Parameter:
1707d763cef2SBarry Smith +  ts     - The TS context
1708d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1709d763cef2SBarry Smith 
1710d763cef2SBarry Smith    Notes:
1711d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1712d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1713d763cef2SBarry Smith    hyphen.
1714d763cef2SBarry Smith 
1715d763cef2SBarry Smith    Contributed by: Matthew Knepley
1716d763cef2SBarry Smith 
1717d763cef2SBarry Smith    Level: advanced
1718d763cef2SBarry Smith 
1719d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1720d763cef2SBarry Smith 
1721d763cef2SBarry Smith .seealso: TSSetFromOptions()
1722d763cef2SBarry Smith 
1723d763cef2SBarry Smith @*/
17241836bdbcSSatish Balay int TSSetOptionsPrefix(TS ts,const char prefix[])
1725d763cef2SBarry Smith {
1726d763cef2SBarry Smith   int ierr;
1727d763cef2SBarry Smith 
1728d763cef2SBarry Smith   PetscFunctionBegin;
1729d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1730d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1731d763cef2SBarry Smith   switch(ts->problem_type) {
1732d763cef2SBarry Smith     case TS_NONLINEAR:
1733d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1734d763cef2SBarry Smith       break;
1735d763cef2SBarry Smith     case TS_LINEAR:
173694b7f48cSBarry Smith       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1737d763cef2SBarry Smith       break;
1738d763cef2SBarry Smith   }
1739d763cef2SBarry Smith   PetscFunctionReturn(0);
1740d763cef2SBarry Smith }
1741d763cef2SBarry Smith 
1742d763cef2SBarry Smith 
17434a2ae208SSatish Balay #undef __FUNCT__
17444a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1745d763cef2SBarry Smith /*@C
1746d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1747d763cef2SBarry Smith    TS options in the database.
1748d763cef2SBarry Smith 
1749d763cef2SBarry Smith    Collective on TS
1750d763cef2SBarry Smith 
1751d763cef2SBarry Smith    Input Parameter:
1752d763cef2SBarry Smith +  ts     - The TS context
1753d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1754d763cef2SBarry Smith 
1755d763cef2SBarry Smith    Notes:
1756d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1757d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1758d763cef2SBarry Smith    hyphen.
1759d763cef2SBarry Smith 
1760d763cef2SBarry Smith    Contributed by: Matthew Knepley
1761d763cef2SBarry Smith 
1762d763cef2SBarry Smith    Level: advanced
1763d763cef2SBarry Smith 
1764d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1765d763cef2SBarry Smith 
1766d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1767d763cef2SBarry Smith 
1768d763cef2SBarry Smith @*/
17691836bdbcSSatish Balay int TSAppendOptionsPrefix(TS ts,const char prefix[])
1770d763cef2SBarry Smith {
1771d763cef2SBarry Smith   int ierr;
1772d763cef2SBarry Smith 
1773d763cef2SBarry Smith   PetscFunctionBegin;
1774d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1775d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1776d763cef2SBarry Smith   switch(ts->problem_type) {
1777d763cef2SBarry Smith     case TS_NONLINEAR:
1778d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1779d763cef2SBarry Smith       break;
1780d763cef2SBarry Smith     case TS_LINEAR:
178194b7f48cSBarry Smith       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1782d763cef2SBarry Smith       break;
1783d763cef2SBarry Smith   }
1784d763cef2SBarry Smith   PetscFunctionReturn(0);
1785d763cef2SBarry Smith }
1786d763cef2SBarry Smith 
17874a2ae208SSatish Balay #undef __FUNCT__
17884a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1789d763cef2SBarry Smith /*@C
1790d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1791d763cef2SBarry Smith    TS options in the database.
1792d763cef2SBarry Smith 
1793d763cef2SBarry Smith    Not Collective
1794d763cef2SBarry Smith 
1795d763cef2SBarry Smith    Input Parameter:
1796d763cef2SBarry Smith .  ts - The TS context
1797d763cef2SBarry Smith 
1798d763cef2SBarry Smith    Output Parameter:
1799d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1800d763cef2SBarry Smith 
1801d763cef2SBarry Smith    Contributed by: Matthew Knepley
1802d763cef2SBarry Smith 
1803d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1804d763cef2SBarry Smith    sufficient length to hold the prefix.
1805d763cef2SBarry Smith 
1806d763cef2SBarry Smith    Level: intermediate
1807d763cef2SBarry Smith 
1808d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1809d763cef2SBarry Smith 
1810d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1811d763cef2SBarry Smith @*/
18121836bdbcSSatish Balay int TSGetOptionsPrefix(TS ts,char *prefix[])
1813d763cef2SBarry Smith {
1814d763cef2SBarry Smith   int ierr;
1815d763cef2SBarry Smith 
1816d763cef2SBarry Smith   PetscFunctionBegin;
1817d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1818d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1819d763cef2SBarry Smith   PetscFunctionReturn(0);
1820d763cef2SBarry Smith }
1821d763cef2SBarry Smith 
18224a2ae208SSatish Balay #undef __FUNCT__
18234a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1824d763cef2SBarry Smith /*@C
1825d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1826d763cef2SBarry Smith 
1827d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1828d763cef2SBarry Smith 
1829d763cef2SBarry Smith    Input Parameter:
1830d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1831d763cef2SBarry Smith 
1832d763cef2SBarry Smith    Output Parameters:
1833d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1834d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1835d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1836d763cef2SBarry Smith 
1837d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1838d763cef2SBarry Smith 
1839d763cef2SBarry Smith    Contributed by: Matthew Knepley
1840d763cef2SBarry Smith 
1841d763cef2SBarry Smith    Level: intermediate
1842d763cef2SBarry Smith 
1843d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1844d763cef2SBarry Smith 
1845d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1846d763cef2SBarry Smith 
1847d763cef2SBarry Smith @*/
1848d763cef2SBarry Smith int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1849d763cef2SBarry Smith {
1850d763cef2SBarry Smith   PetscFunctionBegin;
1851d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1852d763cef2SBarry Smith   if (A)   *A = ts->A;
1853d763cef2SBarry Smith   if (M)   *M = ts->B;
1854d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1855d763cef2SBarry Smith   PetscFunctionReturn(0);
1856d763cef2SBarry Smith }
1857d763cef2SBarry Smith 
18584a2ae208SSatish Balay #undef __FUNCT__
18594a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1860d763cef2SBarry Smith /*@C
1861d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1862d763cef2SBarry Smith 
1863d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1864d763cef2SBarry Smith 
1865d763cef2SBarry Smith    Input Parameter:
1866d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1867d763cef2SBarry Smith 
1868d763cef2SBarry Smith    Output Parameters:
1869d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1870d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1871d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1872d763cef2SBarry Smith 
1873d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1874d763cef2SBarry Smith 
1875d763cef2SBarry Smith    Contributed by: Matthew Knepley
1876d763cef2SBarry Smith 
1877d763cef2SBarry Smith    Level: intermediate
1878d763cef2SBarry Smith 
1879d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1880d763cef2SBarry Smith 
1881d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1882d763cef2SBarry Smith @*/
1883d763cef2SBarry Smith int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1884d763cef2SBarry Smith {
1885d763cef2SBarry Smith   int ierr;
1886d763cef2SBarry Smith 
1887d763cef2SBarry Smith   PetscFunctionBegin;
1888d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1889d763cef2SBarry Smith   PetscFunctionReturn(0);
1890d763cef2SBarry Smith }
1891d763cef2SBarry Smith 
18921713a123SBarry Smith #undef __FUNCT__
18931713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
18941713a123SBarry Smith /*@C
18951713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
18961713a123SBarry Smith    VecView() for the solution at each timestep
18971713a123SBarry Smith 
18981713a123SBarry Smith    Collective on TS
18991713a123SBarry Smith 
19001713a123SBarry Smith    Input Parameters:
19011713a123SBarry Smith +  ts - the TS context
19021713a123SBarry Smith .  step - current time-step
1903142b95e3SSatish Balay .  ptime - current time
19041713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
19051713a123SBarry Smith 
19061713a123SBarry Smith    Level: intermediate
19071713a123SBarry Smith 
19081713a123SBarry Smith .keywords: TS,  vector, monitor, view
19091713a123SBarry Smith 
19101713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
19111713a123SBarry Smith @*/
1912142b95e3SSatish Balay int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
19131713a123SBarry Smith {
19141713a123SBarry Smith   int         ierr;
19151713a123SBarry Smith   PetscViewer viewer = (PetscViewer) dummy;
19161713a123SBarry Smith 
19171713a123SBarry Smith   PetscFunctionBegin;
19181713a123SBarry Smith   if (!viewer) {
19191713a123SBarry Smith     MPI_Comm comm;
19201713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
19211713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
19221713a123SBarry Smith   }
19231713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
19241713a123SBarry Smith   PetscFunctionReturn(0);
19251713a123SBarry Smith }
19261713a123SBarry Smith 
19271713a123SBarry Smith 
19281713a123SBarry Smith 
1929