xref: /petsc/src/ts/interface/ts.c (revision 3c94ec114b93ba0cfbf7331735f0e93b9db506db)
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;
262fc52814SBarry 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;
804482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
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;
2374482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2384482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
239c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
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 */
2514482741eSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE,4);
2524482741eSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE,5);
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;
2774482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2784482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
2794482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
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 
3384482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
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;
3964482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
3974482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
3984482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
399c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
400c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
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.
418*3c94ec11SBarry Smith    Use TSSetRHSMatrix() for linear problems.
419d763cef2SBarry Smith 
420d763cef2SBarry Smith    Collective on TS
421d763cef2SBarry Smith 
422d763cef2SBarry Smith    Input Parameters:
423d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
424d763cef2SBarry Smith .  A   - Jacobian matrix
425d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
426d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
427d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
428d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
429d763cef2SBarry Smith 
430d763cef2SBarry Smith    Calling sequence of func:
43187828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
432d763cef2SBarry Smith 
433d763cef2SBarry Smith +  t - current timestep
434d763cef2SBarry Smith .  u - input vector
435d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
436d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
437d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
43894b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
439d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
440d763cef2SBarry Smith 
441d763cef2SBarry Smith    Notes:
44294b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
443d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
444d763cef2SBarry Smith 
445d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
446d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
447d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
448d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
449d763cef2SBarry Smith    throughout the global iterations.
450d763cef2SBarry Smith 
451d763cef2SBarry Smith    Level: beginner
452d763cef2SBarry Smith 
453d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
454d763cef2SBarry Smith 
455d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
456*3c94ec11SBarry Smith           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetRHSMatrix()
457d763cef2SBarry Smith 
458d763cef2SBarry Smith @*/
45987828ca2SBarry Smith int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
460d763cef2SBarry Smith {
461d763cef2SBarry Smith   PetscFunctionBegin;
4624482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4634482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
4644482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
465c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
466c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
467d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
46829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
469d763cef2SBarry Smith   }
470d763cef2SBarry Smith 
471000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
472d763cef2SBarry Smith   ts->jacP             = ctx;
473d763cef2SBarry Smith   ts->A                = A;
474d763cef2SBarry Smith   ts->B                = B;
475d763cef2SBarry Smith   PetscFunctionReturn(0);
476d763cef2SBarry Smith }
477d763cef2SBarry Smith 
4784a2ae208SSatish Balay #undef __FUNCT__
4794a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
480d763cef2SBarry Smith /*
481d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
482d763cef2SBarry Smith 
483d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
484d763cef2SBarry Smith    this routine applies the matrix.
485d763cef2SBarry Smith */
48687828ca2SBarry Smith int TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
487d763cef2SBarry Smith {
488d763cef2SBarry Smith   int ierr;
489d763cef2SBarry Smith 
490d763cef2SBarry Smith   PetscFunctionBegin;
4914482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4924482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
493c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,x,3);
494d763cef2SBarry Smith 
495000e7ae3SMatthew Knepley   if (ts->ops->rhsbc) {
496d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
497000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
498d763cef2SBarry Smith     PetscStackPop;
499d763cef2SBarry Smith     PetscFunctionReturn(0);
500d763cef2SBarry Smith   }
501d763cef2SBarry Smith 
502d763cef2SBarry Smith   PetscFunctionReturn(0);
503d763cef2SBarry Smith }
504d763cef2SBarry Smith 
5054a2ae208SSatish Balay #undef __FUNCT__
5064a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
507d763cef2SBarry Smith /*@C
508d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
509d763cef2SBarry Smith     boundary conditions for the function F.
510d763cef2SBarry Smith 
511d763cef2SBarry Smith     Collective on TS
512d763cef2SBarry Smith 
513d763cef2SBarry Smith     Input Parameters:
514d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
515d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
516d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
517d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
518d763cef2SBarry Smith 
519d763cef2SBarry Smith     Calling sequence of func:
52087828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
521d763cef2SBarry Smith 
522d763cef2SBarry Smith +   t - current timestep
523d763cef2SBarry Smith .   F - function vector
524d763cef2SBarry Smith -   ctx - [optional] user-defined function context
525d763cef2SBarry Smith 
526d763cef2SBarry Smith     Level: intermediate
527d763cef2SBarry Smith 
528d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
529d763cef2SBarry Smith @*/
53087828ca2SBarry Smith int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,PetscReal,Vec,void*),void *ctx)
531d763cef2SBarry Smith {
532d763cef2SBarry Smith   PetscFunctionBegin;
533d763cef2SBarry Smith 
5344482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
535d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
53629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
537d763cef2SBarry Smith   }
538000e7ae3SMatthew Knepley   ts->ops->rhsbc = f;
539d763cef2SBarry Smith   ts->bcP        = ctx;
540d763cef2SBarry Smith   PetscFunctionReturn(0);
541d763cef2SBarry Smith }
542d763cef2SBarry Smith 
5434a2ae208SSatish Balay #undef __FUNCT__
5444a2ae208SSatish Balay #define __FUNCT__ "TSView"
5457e2c5f70SBarry Smith /*@C
546d763cef2SBarry Smith     TSView - Prints the TS data structure.
547d763cef2SBarry Smith 
5484c49b128SBarry Smith     Collective on TS
549d763cef2SBarry Smith 
550d763cef2SBarry Smith     Input Parameters:
551d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
552d763cef2SBarry Smith -   viewer - visualization context
553d763cef2SBarry Smith 
554d763cef2SBarry Smith     Options Database Key:
555d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
556d763cef2SBarry Smith 
557d763cef2SBarry Smith     Notes:
558d763cef2SBarry Smith     The available visualization contexts include
559b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
560b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
561d763cef2SBarry Smith          output where only the first processor opens
562d763cef2SBarry Smith          the file.  All other processors send their
563d763cef2SBarry Smith          data to the first processor to print.
564d763cef2SBarry Smith 
565d763cef2SBarry Smith     The user can open an alternative visualization context with
566b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
567d763cef2SBarry Smith 
568d763cef2SBarry Smith     Level: beginner
569d763cef2SBarry Smith 
570d763cef2SBarry Smith .keywords: TS, timestep, view
571d763cef2SBarry Smith 
572b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
573d763cef2SBarry Smith @*/
574b0a32e0cSBarry Smith int TSView(TS ts,PetscViewer viewer)
575d763cef2SBarry Smith {
576d763cef2SBarry Smith   int        ierr;
577454a90a3SBarry Smith   char       *type;
5786831982aSBarry Smith   PetscTruth isascii,isstring;
579d763cef2SBarry Smith 
580d763cef2SBarry Smith   PetscFunctionBegin;
5814482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
582b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
5834482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
584c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
585fd16b177SBarry Smith 
586b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
587b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
5880f5bd95cSBarry Smith   if (isascii) {
589b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
590454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
591454a90a3SBarry Smith     if (type) {
592b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
593184914b5SBarry Smith     } else {
594b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
595184914b5SBarry Smith     }
596000e7ae3SMatthew Knepley     if (ts->ops->view) {
597b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
598000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
599b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
600d763cef2SBarry Smith     }
601b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
602b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
603d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
604b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
605d763cef2SBarry Smith     }
606b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
6070f5bd95cSBarry Smith   } else if (isstring) {
608454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
609b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
610d763cef2SBarry Smith   }
611b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
61294b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
613d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
614b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
615d763cef2SBarry Smith   PetscFunctionReturn(0);
616d763cef2SBarry Smith }
617d763cef2SBarry Smith 
618d763cef2SBarry Smith 
6194a2ae208SSatish Balay #undef __FUNCT__
6204a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
621d763cef2SBarry Smith /*@C
622d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
623d763cef2SBarry Smith    the timesteppers.
624d763cef2SBarry Smith 
625d763cef2SBarry Smith    Collective on TS
626d763cef2SBarry Smith 
627d763cef2SBarry Smith    Input Parameters:
628d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
629d763cef2SBarry Smith -  usrP - optional user context
630d763cef2SBarry Smith 
631d763cef2SBarry Smith    Level: intermediate
632d763cef2SBarry Smith 
633d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
634d763cef2SBarry Smith 
635d763cef2SBarry Smith .seealso: TSGetApplicationContext()
636d763cef2SBarry Smith @*/
637d763cef2SBarry Smith int TSSetApplicationContext(TS ts,void *usrP)
638d763cef2SBarry Smith {
639d763cef2SBarry Smith   PetscFunctionBegin;
6404482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
641d763cef2SBarry Smith   ts->user = usrP;
642d763cef2SBarry Smith   PetscFunctionReturn(0);
643d763cef2SBarry Smith }
644d763cef2SBarry Smith 
6454a2ae208SSatish Balay #undef __FUNCT__
6464a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
647d763cef2SBarry Smith /*@C
648d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
649d763cef2SBarry Smith     timestepper.
650d763cef2SBarry Smith 
651d763cef2SBarry Smith     Not Collective
652d763cef2SBarry Smith 
653d763cef2SBarry Smith     Input Parameter:
654d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
655d763cef2SBarry Smith 
656d763cef2SBarry Smith     Output Parameter:
657d763cef2SBarry Smith .   usrP - user context
658d763cef2SBarry Smith 
659d763cef2SBarry Smith     Level: intermediate
660d763cef2SBarry Smith 
661d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
662d763cef2SBarry Smith 
663d763cef2SBarry Smith .seealso: TSSetApplicationContext()
664d763cef2SBarry Smith @*/
665d763cef2SBarry Smith int TSGetApplicationContext(TS ts,void **usrP)
666d763cef2SBarry Smith {
667d763cef2SBarry Smith   PetscFunctionBegin;
6684482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
669d763cef2SBarry Smith   *usrP = ts->user;
670d763cef2SBarry Smith   PetscFunctionReturn(0);
671d763cef2SBarry Smith }
672d763cef2SBarry Smith 
6734a2ae208SSatish Balay #undef __FUNCT__
6744a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
675d763cef2SBarry Smith /*@
676d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
677d763cef2SBarry Smith 
678d763cef2SBarry Smith    Not Collective
679d763cef2SBarry Smith 
680d763cef2SBarry Smith    Input Parameter:
681d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
682d763cef2SBarry Smith 
683d763cef2SBarry Smith    Output Parameter:
684d763cef2SBarry Smith .  iter - number steps so far
685d763cef2SBarry Smith 
686d763cef2SBarry Smith    Level: intermediate
687d763cef2SBarry Smith 
688d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
689d763cef2SBarry Smith @*/
690d763cef2SBarry Smith int TSGetTimeStepNumber(TS ts,int* iter)
691d763cef2SBarry Smith {
692d763cef2SBarry Smith   PetscFunctionBegin;
6934482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
6944482741eSBarry Smith   PetscValidIntPointer(iter,2);
695d763cef2SBarry Smith   *iter = ts->steps;
696d763cef2SBarry Smith   PetscFunctionReturn(0);
697d763cef2SBarry Smith }
698d763cef2SBarry Smith 
6994a2ae208SSatish Balay #undef __FUNCT__
7004a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
701d763cef2SBarry Smith /*@
702d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
703d763cef2SBarry Smith    as well as the initial time.
704d763cef2SBarry Smith 
705d763cef2SBarry Smith    Collective on TS
706d763cef2SBarry Smith 
707d763cef2SBarry Smith    Input Parameters:
708d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
709d763cef2SBarry Smith .  initial_time - the initial time
710d763cef2SBarry Smith -  time_step - the size of the timestep
711d763cef2SBarry Smith 
712d763cef2SBarry Smith    Level: intermediate
713d763cef2SBarry Smith 
714d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
715d763cef2SBarry Smith 
716d763cef2SBarry Smith .keywords: TS, set, initial, timestep
717d763cef2SBarry Smith @*/
71887828ca2SBarry Smith int TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
719d763cef2SBarry Smith {
720d763cef2SBarry Smith   PetscFunctionBegin;
7214482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
722d763cef2SBarry Smith   ts->time_step         = time_step;
723d763cef2SBarry Smith   ts->initial_time_step = time_step;
724d763cef2SBarry Smith   ts->ptime             = initial_time;
725d763cef2SBarry Smith   PetscFunctionReturn(0);
726d763cef2SBarry Smith }
727d763cef2SBarry Smith 
7284a2ae208SSatish Balay #undef __FUNCT__
7294a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
730d763cef2SBarry Smith /*@
731d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
732d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
733d763cef2SBarry Smith 
734d763cef2SBarry Smith    Collective on TS
735d763cef2SBarry Smith 
736d763cef2SBarry Smith    Input Parameters:
737d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
738d763cef2SBarry Smith -  time_step - the size of the timestep
739d763cef2SBarry Smith 
740d763cef2SBarry Smith    Level: intermediate
741d763cef2SBarry Smith 
742d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
743d763cef2SBarry Smith 
744d763cef2SBarry Smith .keywords: TS, set, timestep
745d763cef2SBarry Smith @*/
74687828ca2SBarry Smith int TSSetTimeStep(TS ts,PetscReal time_step)
747d763cef2SBarry Smith {
748d763cef2SBarry Smith   PetscFunctionBegin;
7494482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
750d763cef2SBarry Smith   ts->time_step = time_step;
751d763cef2SBarry Smith   PetscFunctionReturn(0);
752d763cef2SBarry Smith }
753d763cef2SBarry Smith 
7544a2ae208SSatish Balay #undef __FUNCT__
7554a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
756d763cef2SBarry Smith /*@
757d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
758d763cef2SBarry Smith 
759d763cef2SBarry Smith    Not Collective
760d763cef2SBarry Smith 
761d763cef2SBarry Smith    Input Parameter:
762d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
763d763cef2SBarry Smith 
764d763cef2SBarry Smith    Output Parameter:
765d763cef2SBarry Smith .  dt - the current timestep size
766d763cef2SBarry Smith 
767d763cef2SBarry Smith    Level: intermediate
768d763cef2SBarry Smith 
769d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
770d763cef2SBarry Smith 
771d763cef2SBarry Smith .keywords: TS, get, timestep
772d763cef2SBarry Smith @*/
77387828ca2SBarry Smith int TSGetTimeStep(TS ts,PetscReal* dt)
774d763cef2SBarry Smith {
775d763cef2SBarry Smith   PetscFunctionBegin;
7764482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
7774482741eSBarry Smith   PetscValidDoublePointer(dt,2);
778d763cef2SBarry Smith   *dt = ts->time_step;
779d763cef2SBarry Smith   PetscFunctionReturn(0);
780d763cef2SBarry Smith }
781d763cef2SBarry Smith 
7824a2ae208SSatish Balay #undef __FUNCT__
7834a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
784d763cef2SBarry Smith /*@C
785d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
786d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
787d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
788d763cef2SBarry Smith    the solution at the next timestep has been calculated.
789d763cef2SBarry Smith 
790d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
791d763cef2SBarry Smith 
792d763cef2SBarry Smith    Input Parameter:
793d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
794d763cef2SBarry Smith 
795d763cef2SBarry Smith    Output Parameter:
796d763cef2SBarry Smith .  v - the vector containing the solution
797d763cef2SBarry Smith 
798d763cef2SBarry Smith    Level: intermediate
799d763cef2SBarry Smith 
800d763cef2SBarry Smith .seealso: TSGetTimeStep()
801d763cef2SBarry Smith 
802d763cef2SBarry Smith .keywords: TS, timestep, get, solution
803d763cef2SBarry Smith @*/
804d763cef2SBarry Smith int TSGetSolution(TS ts,Vec *v)
805d763cef2SBarry Smith {
806d763cef2SBarry Smith   PetscFunctionBegin;
8074482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
8084482741eSBarry Smith   PetscValidPointer(v,2);
809d763cef2SBarry Smith   *v = ts->vec_sol_always;
810d763cef2SBarry Smith   PetscFunctionReturn(0);
811d763cef2SBarry Smith }
812d763cef2SBarry Smith 
813bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
8144a2ae208SSatish Balay #undef __FUNCT__
815bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
816d763cef2SBarry Smith /*@C
817bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
818d763cef2SBarry Smith 
819bdad233fSMatthew Knepley   Not collective
820d763cef2SBarry Smith 
821bdad233fSMatthew Knepley   Input Parameters:
822bdad233fSMatthew Knepley + ts   - The TS
823bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
824d763cef2SBarry Smith .vb
825d763cef2SBarry Smith          U_t = A U
826d763cef2SBarry Smith          U_t = A(t) U
827d763cef2SBarry Smith          U_t = F(t,U)
828d763cef2SBarry Smith .ve
829d763cef2SBarry Smith 
830d763cef2SBarry Smith    Level: beginner
831d763cef2SBarry Smith 
832bdad233fSMatthew Knepley .keywords: TS, problem type
833bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
834d763cef2SBarry Smith @*/
835bdad233fSMatthew Knepley int TSSetProblemType(TS ts, TSProblemType type) {
836d763cef2SBarry Smith   PetscFunctionBegin;
8374482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
838bdad233fSMatthew Knepley   ts->problem_type = type;
839d763cef2SBarry Smith   PetscFunctionReturn(0);
840d763cef2SBarry Smith }
841d763cef2SBarry Smith 
842bdad233fSMatthew Knepley #undef __FUNCT__
843bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
844bdad233fSMatthew Knepley /*@C
845bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
846bdad233fSMatthew Knepley 
847bdad233fSMatthew Knepley   Not collective
848bdad233fSMatthew Knepley 
849bdad233fSMatthew Knepley   Input Parameter:
850bdad233fSMatthew Knepley . ts   - The TS
851bdad233fSMatthew Knepley 
852bdad233fSMatthew Knepley   Output Parameter:
853bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
854bdad233fSMatthew Knepley .vb
855bdad233fSMatthew Knepley          U_t = A U
856bdad233fSMatthew Knepley          U_t = A(t) U
857bdad233fSMatthew Knepley          U_t = F(t,U)
858bdad233fSMatthew Knepley .ve
859bdad233fSMatthew Knepley 
860bdad233fSMatthew Knepley    Level: beginner
861bdad233fSMatthew Knepley 
862bdad233fSMatthew Knepley .keywords: TS, problem type
863bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
864bdad233fSMatthew Knepley @*/
865bdad233fSMatthew Knepley int TSGetProblemType(TS ts, TSProblemType *type) {
866bdad233fSMatthew Knepley   PetscFunctionBegin;
8674482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
8684482741eSBarry Smith   PetscValidIntPointer(type,2);
869bdad233fSMatthew Knepley   *type = ts->problem_type;
870bdad233fSMatthew Knepley   PetscFunctionReturn(0);
871bdad233fSMatthew Knepley }
872d763cef2SBarry Smith 
8734a2ae208SSatish Balay #undef __FUNCT__
8744a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
875d763cef2SBarry Smith /*@
876d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
877d763cef2SBarry Smith    of a timestepper.
878d763cef2SBarry Smith 
879d763cef2SBarry Smith    Collective on TS
880d763cef2SBarry Smith 
881d763cef2SBarry Smith    Input Parameter:
882d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
883d763cef2SBarry Smith 
884d763cef2SBarry Smith    Notes:
885d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
886d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
887d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
888d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
889d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
890d763cef2SBarry Smith 
891d763cef2SBarry Smith    Level: advanced
892d763cef2SBarry Smith 
893d763cef2SBarry Smith .keywords: TS, timestep, setup
894d763cef2SBarry Smith 
895d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
896d763cef2SBarry Smith @*/
897d763cef2SBarry Smith int TSSetUp(TS ts)
898d763cef2SBarry Smith {
899d763cef2SBarry Smith   int ierr;
900d763cef2SBarry Smith 
901d763cef2SBarry Smith   PetscFunctionBegin;
9024482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
90329bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
904d763cef2SBarry Smith   if (!ts->type_name) {
905d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
906d763cef2SBarry Smith   }
907000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
908d763cef2SBarry Smith   ts->setupcalled = 1;
909d763cef2SBarry Smith   PetscFunctionReturn(0);
910d763cef2SBarry Smith }
911d763cef2SBarry Smith 
9124a2ae208SSatish Balay #undef __FUNCT__
9134a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
914d763cef2SBarry Smith /*@C
915d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
916d763cef2SBarry Smith    with TSCreate().
917d763cef2SBarry Smith 
918d763cef2SBarry Smith    Collective on TS
919d763cef2SBarry Smith 
920d763cef2SBarry Smith    Input Parameter:
921d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
922d763cef2SBarry Smith 
923d763cef2SBarry Smith    Level: beginner
924d763cef2SBarry Smith 
925d763cef2SBarry Smith .keywords: TS, timestepper, destroy
926d763cef2SBarry Smith 
927d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
928d763cef2SBarry Smith @*/
929d763cef2SBarry Smith int TSDestroy(TS ts)
930d763cef2SBarry Smith {
931329f5518SBarry Smith   int ierr,i;
932d763cef2SBarry Smith 
933d763cef2SBarry Smith   PetscFunctionBegin;
9344482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
935d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
936d763cef2SBarry Smith 
937be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
9380f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
939be0abb6dSBarry Smith 
94094b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
941d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
942000e7ae3SMatthew Knepley   ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);
943329f5518SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
944329f5518SBarry Smith     if (ts->mdestroy[i]) {
945329f5518SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
946329f5518SBarry Smith     }
947329f5518SBarry Smith   }
948b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)ts);
949d763cef2SBarry Smith   PetscHeaderDestroy((PetscObject)ts);
950d763cef2SBarry Smith   PetscFunctionReturn(0);
951d763cef2SBarry Smith }
952d763cef2SBarry Smith 
9534a2ae208SSatish Balay #undef __FUNCT__
9544a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
955d763cef2SBarry Smith /*@C
956d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
957d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
958d763cef2SBarry Smith 
959d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
960d763cef2SBarry Smith 
961d763cef2SBarry Smith    Input Parameter:
962d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
963d763cef2SBarry Smith 
964d763cef2SBarry Smith    Output Parameter:
965d763cef2SBarry Smith .  snes - the nonlinear solver context
966d763cef2SBarry Smith 
967d763cef2SBarry Smith    Notes:
968d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
969d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
97094b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
971d763cef2SBarry Smith 
972d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
973d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
974d763cef2SBarry Smith 
975d763cef2SBarry Smith    Level: beginner
976d763cef2SBarry Smith 
977d763cef2SBarry Smith .keywords: timestep, get, SNES
978d763cef2SBarry Smith @*/
979d763cef2SBarry Smith int TSGetSNES(TS ts,SNES *snes)
980d763cef2SBarry Smith {
981d763cef2SBarry Smith   PetscFunctionBegin;
9824482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
9834482741eSBarry Smith   PetscValidPointer(snes,2);
98494b7f48cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
985d763cef2SBarry Smith   *snes = ts->snes;
986d763cef2SBarry Smith   PetscFunctionReturn(0);
987d763cef2SBarry Smith }
988d763cef2SBarry Smith 
9894a2ae208SSatish Balay #undef __FUNCT__
99094b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
991d763cef2SBarry Smith /*@C
99294b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
993d763cef2SBarry Smith    a TS (timestepper) context.
994d763cef2SBarry Smith 
99594b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
996d763cef2SBarry Smith 
997d763cef2SBarry Smith    Input Parameter:
998d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
999d763cef2SBarry Smith 
1000d763cef2SBarry Smith    Output Parameter:
100194b7f48cSBarry Smith .  ksp - the nonlinear solver context
1002d763cef2SBarry Smith 
1003d763cef2SBarry Smith    Notes:
100494b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1005d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1006d763cef2SBarry Smith    KSP and PC contexts as well.
1007d763cef2SBarry Smith 
100894b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
100994b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1010d763cef2SBarry Smith 
1011d763cef2SBarry Smith    Level: beginner
1012d763cef2SBarry Smith 
101394b7f48cSBarry Smith .keywords: timestep, get, KSP
1014d763cef2SBarry Smith @*/
101594b7f48cSBarry Smith int TSGetKSP(TS ts,KSP *ksp)
1016d763cef2SBarry Smith {
1017d763cef2SBarry Smith   PetscFunctionBegin;
10184482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
10194482741eSBarry Smith   PetscValidPointer(ksp,2);
102029bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
102194b7f48cSBarry Smith   *ksp = ts->ksp;
1022d763cef2SBarry Smith   PetscFunctionReturn(0);
1023d763cef2SBarry Smith }
1024d763cef2SBarry Smith 
1025d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1026d763cef2SBarry Smith 
10274a2ae208SSatish Balay #undef __FUNCT__
1028adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1029adb62b0dSMatthew Knepley /*@
1030adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1031adb62b0dSMatthew Knepley    maximum time for iteration.
1032adb62b0dSMatthew Knepley 
1033adb62b0dSMatthew Knepley    Collective on TS
1034adb62b0dSMatthew Knepley 
1035adb62b0dSMatthew Knepley    Input Parameters:
1036adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1037adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1038adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1039adb62b0dSMatthew Knepley 
1040adb62b0dSMatthew Knepley    Level: intermediate
1041adb62b0dSMatthew Knepley 
1042adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1043adb62b0dSMatthew Knepley @*/
1044adb62b0dSMatthew Knepley int TSGetDuration(TS ts, int *maxsteps, PetscReal *maxtime)
1045adb62b0dSMatthew Knepley {
1046adb62b0dSMatthew Knepley   PetscFunctionBegin;
10474482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1048adb62b0dSMatthew Knepley   if (maxsteps != PETSC_NULL) {
10494482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1050adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1051adb62b0dSMatthew Knepley   }
1052adb62b0dSMatthew Knepley   if (maxtime  != PETSC_NULL) {
10534482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1054adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1055adb62b0dSMatthew Knepley   }
1056adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1057adb62b0dSMatthew Knepley }
1058adb62b0dSMatthew Knepley 
1059adb62b0dSMatthew Knepley #undef __FUNCT__
10604a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1061d763cef2SBarry Smith /*@
1062d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1063d763cef2SBarry Smith    maximum time for iteration.
1064d763cef2SBarry Smith 
1065d763cef2SBarry Smith    Collective on TS
1066d763cef2SBarry Smith 
1067d763cef2SBarry Smith    Input Parameters:
1068d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1069d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1070d763cef2SBarry Smith -  maxtime - final time to iterate to
1071d763cef2SBarry Smith 
1072d763cef2SBarry Smith    Options Database Keys:
1073d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1074d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1075d763cef2SBarry Smith 
1076d763cef2SBarry Smith    Notes:
1077d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1078d763cef2SBarry Smith 
1079d763cef2SBarry Smith    Level: intermediate
1080d763cef2SBarry Smith 
1081d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1082d763cef2SBarry Smith @*/
108387828ca2SBarry Smith int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime)
1084d763cef2SBarry Smith {
1085d763cef2SBarry Smith   PetscFunctionBegin;
10864482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1087d763cef2SBarry Smith   ts->max_steps = maxsteps;
1088d763cef2SBarry Smith   ts->max_time  = maxtime;
1089d763cef2SBarry Smith   PetscFunctionReturn(0);
1090d763cef2SBarry Smith }
1091d763cef2SBarry Smith 
10924a2ae208SSatish Balay #undef __FUNCT__
10934a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1094d763cef2SBarry Smith /*@
1095d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1096d763cef2SBarry Smith    for use by the TS routines.
1097d763cef2SBarry Smith 
1098d763cef2SBarry Smith    Collective on TS and Vec
1099d763cef2SBarry Smith 
1100d763cef2SBarry Smith    Input Parameters:
1101d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1102d763cef2SBarry Smith -  x - the solution vector
1103d763cef2SBarry Smith 
1104d763cef2SBarry Smith    Level: beginner
1105d763cef2SBarry Smith 
1106d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1107d763cef2SBarry Smith @*/
1108d763cef2SBarry Smith int TSSetSolution(TS ts,Vec x)
1109d763cef2SBarry Smith {
1110d763cef2SBarry Smith   PetscFunctionBegin;
11114482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
11124482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1113d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1114d763cef2SBarry Smith   PetscFunctionReturn(0);
1115d763cef2SBarry Smith }
1116d763cef2SBarry Smith 
1117e74ef692SMatthew Knepley #undef __FUNCT__
1118e74ef692SMatthew Knepley #define __FUNCT__ "TSSetRhsBC"
1119ac226902SBarry Smith /*@C
1120000e7ae3SMatthew Knepley   TSSetRhsBC - Sets the function which applies boundary conditions
1121000e7ae3SMatthew Knepley   to the Rhs of each system.
1122000e7ae3SMatthew Knepley 
1123000e7ae3SMatthew Knepley   Collective on TS
1124000e7ae3SMatthew Knepley 
1125000e7ae3SMatthew Knepley   Input Parameters:
1126000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1127000e7ae3SMatthew Knepley - func - The function
1128000e7ae3SMatthew Knepley 
1129000e7ae3SMatthew Knepley   Calling sequence of func:
1130000e7ae3SMatthew Knepley . func (TS ts, Vec rhs, void *ctx);
1131000e7ae3SMatthew Knepley 
1132000e7ae3SMatthew Knepley + rhs - The current rhs vector
1133000e7ae3SMatthew Knepley - ctx - The user-context
1134000e7ae3SMatthew Knepley 
1135000e7ae3SMatthew Knepley   Level: intermediate
1136000e7ae3SMatthew Knepley 
1137000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1138000e7ae3SMatthew Knepley @*/
1139000e7ae3SMatthew Knepley int TSSetRhsBC(TS ts, int (*func)(TS, Vec, void *))
1140000e7ae3SMatthew Knepley {
1141000e7ae3SMatthew Knepley   PetscFunctionBegin;
11424482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1143000e7ae3SMatthew Knepley   ts->ops->applyrhsbc = func;
1144000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1145000e7ae3SMatthew Knepley }
1146000e7ae3SMatthew Knepley 
1147e74ef692SMatthew Knepley #undef __FUNCT__
1148e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultRhsBC"
1149000e7ae3SMatthew Knepley /*@
1150000e7ae3SMatthew Knepley   TSDefaultRhsBC - The default boundary condition function which does nothing.
1151000e7ae3SMatthew Knepley 
1152000e7ae3SMatthew Knepley   Collective on TS
1153000e7ae3SMatthew Knepley 
1154000e7ae3SMatthew Knepley   Input Parameters:
1155000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1156000e7ae3SMatthew Knepley . rhs - The Rhs
1157000e7ae3SMatthew Knepley - ctx - The user-context
1158000e7ae3SMatthew Knepley 
1159000e7ae3SMatthew Knepley   Level: developer
1160000e7ae3SMatthew Knepley 
1161000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1162000e7ae3SMatthew Knepley @*/
1163000e7ae3SMatthew Knepley int TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1164000e7ae3SMatthew Knepley {
1165000e7ae3SMatthew Knepley   PetscFunctionBegin;
1166000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1167000e7ae3SMatthew Knepley }
1168000e7ae3SMatthew Knepley 
1169e74ef692SMatthew Knepley #undef __FUNCT__
1170e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSystemMatrixBC"
1171ac226902SBarry Smith /*@C
1172000e7ae3SMatthew Knepley   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1173000e7ae3SMatthew Knepley   to the system matrix and preconditioner of each system.
1174000e7ae3SMatthew Knepley 
1175000e7ae3SMatthew Knepley   Collective on TS
1176000e7ae3SMatthew Knepley 
1177000e7ae3SMatthew Knepley   Input Parameters:
1178000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1179000e7ae3SMatthew Knepley - func - The function
1180000e7ae3SMatthew Knepley 
1181000e7ae3SMatthew Knepley   Calling sequence of func:
1182000e7ae3SMatthew Knepley . func (TS ts, Mat A, Mat B, void *ctx);
1183000e7ae3SMatthew Knepley 
1184000e7ae3SMatthew Knepley + A   - The current system matrix
1185000e7ae3SMatthew Knepley . B   - The current preconditioner
1186000e7ae3SMatthew Knepley - ctx - The user-context
1187000e7ae3SMatthew Knepley 
1188000e7ae3SMatthew Knepley   Level: intermediate
1189000e7ae3SMatthew Knepley 
1190000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1191000e7ae3SMatthew Knepley @*/
1192000e7ae3SMatthew Knepley int TSSetSystemMatrixBC(TS ts, int (*func)(TS, Mat, Mat, void *))
1193000e7ae3SMatthew Knepley {
1194000e7ae3SMatthew Knepley   PetscFunctionBegin;
11954482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1196000e7ae3SMatthew Knepley   ts->ops->applymatrixbc = func;
1197000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1198000e7ae3SMatthew Knepley }
1199000e7ae3SMatthew Knepley 
1200e74ef692SMatthew Knepley #undef __FUNCT__
1201e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSystemMatrixBC"
1202000e7ae3SMatthew Knepley /*@
1203000e7ae3SMatthew Knepley   TSDefaultSystemMatrixBC - The default boundary condition function which
1204000e7ae3SMatthew Knepley   does nothing.
1205000e7ae3SMatthew Knepley 
1206000e7ae3SMatthew Knepley   Collective on TS
1207000e7ae3SMatthew Knepley 
1208000e7ae3SMatthew Knepley   Input Parameters:
1209000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1210000e7ae3SMatthew Knepley . A   - The system matrix
1211000e7ae3SMatthew Knepley . B   - The preconditioner
1212000e7ae3SMatthew Knepley - ctx - The user-context
1213000e7ae3SMatthew Knepley 
1214000e7ae3SMatthew Knepley   Level: developer
1215000e7ae3SMatthew Knepley 
1216000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1217000e7ae3SMatthew Knepley @*/
1218000e7ae3SMatthew Knepley int TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1219000e7ae3SMatthew Knepley {
1220000e7ae3SMatthew Knepley   PetscFunctionBegin;
1221000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1222000e7ae3SMatthew Knepley }
1223000e7ae3SMatthew Knepley 
1224e74ef692SMatthew Knepley #undef __FUNCT__
1225e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSolutionBC"
1226ac226902SBarry Smith /*@C
1227000e7ae3SMatthew Knepley   TSSetSolutionBC - Sets the function which applies boundary conditions
1228000e7ae3SMatthew Knepley   to the solution of each system. This is necessary in nonlinear systems
1229000e7ae3SMatthew Knepley   which time dependent boundary conditions.
1230000e7ae3SMatthew Knepley 
1231000e7ae3SMatthew Knepley   Collective on TS
1232000e7ae3SMatthew Knepley 
1233000e7ae3SMatthew Knepley   Input Parameters:
1234000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1235000e7ae3SMatthew Knepley - func - The function
1236000e7ae3SMatthew Knepley 
1237000e7ae3SMatthew Knepley   Calling sequence of func:
1238000e7ae3SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
1239000e7ae3SMatthew Knepley 
1240000e7ae3SMatthew Knepley + sol - The current solution vector
1241000e7ae3SMatthew Knepley - ctx - The user-context
1242000e7ae3SMatthew Knepley 
1243000e7ae3SMatthew Knepley   Level: intermediate
1244000e7ae3SMatthew Knepley 
1245000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1246000e7ae3SMatthew Knepley @*/
1247000e7ae3SMatthew Knepley int TSSetSolutionBC(TS ts, int (*func)(TS, Vec, void *))
1248000e7ae3SMatthew Knepley {
1249000e7ae3SMatthew Knepley   PetscFunctionBegin;
12504482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1251000e7ae3SMatthew Knepley   ts->ops->applysolbc = func;
1252000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1253000e7ae3SMatthew Knepley }
1254000e7ae3SMatthew Knepley 
1255e74ef692SMatthew Knepley #undef __FUNCT__
1256e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSolutionBC"
1257000e7ae3SMatthew Knepley /*@
1258000e7ae3SMatthew Knepley   TSDefaultSolutionBC - The default boundary condition function which
1259000e7ae3SMatthew Knepley   does nothing.
1260000e7ae3SMatthew Knepley 
1261000e7ae3SMatthew Knepley   Collective on TS
1262000e7ae3SMatthew Knepley 
1263000e7ae3SMatthew Knepley   Input Parameters:
1264000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1265000e7ae3SMatthew Knepley . sol - The solution
1266000e7ae3SMatthew Knepley - ctx - The user-context
1267000e7ae3SMatthew Knepley 
1268000e7ae3SMatthew Knepley   Level: developer
1269000e7ae3SMatthew Knepley 
1270000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1271000e7ae3SMatthew Knepley @*/
1272000e7ae3SMatthew Knepley int TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1273000e7ae3SMatthew Knepley {
1274000e7ae3SMatthew Knepley   PetscFunctionBegin;
1275000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1276000e7ae3SMatthew Knepley }
1277000e7ae3SMatthew Knepley 
1278e74ef692SMatthew Knepley #undef __FUNCT__
1279e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1280ac226902SBarry Smith /*@C
1281000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1282000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1283000e7ae3SMatthew Knepley 
1284000e7ae3SMatthew Knepley   Collective on TS
1285000e7ae3SMatthew Knepley 
1286000e7ae3SMatthew Knepley   Input Parameters:
1287000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1288000e7ae3SMatthew Knepley - func - The function
1289000e7ae3SMatthew Knepley 
1290000e7ae3SMatthew Knepley   Calling sequence of func:
1291000e7ae3SMatthew Knepley . func (TS ts);
1292000e7ae3SMatthew Knepley 
1293000e7ae3SMatthew Knepley   Level: intermediate
1294000e7ae3SMatthew Knepley 
1295000e7ae3SMatthew Knepley .keywords: TS, timestep
1296000e7ae3SMatthew Knepley @*/
1297000e7ae3SMatthew Knepley int TSSetPreStep(TS ts, int (*func)(TS))
1298000e7ae3SMatthew Knepley {
1299000e7ae3SMatthew Knepley   PetscFunctionBegin;
13004482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1301000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1302000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1303000e7ae3SMatthew Knepley }
1304000e7ae3SMatthew Knepley 
1305e74ef692SMatthew Knepley #undef __FUNCT__
1306e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1307000e7ae3SMatthew Knepley /*@
1308000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1309000e7ae3SMatthew Knepley 
1310000e7ae3SMatthew Knepley   Collective on TS
1311000e7ae3SMatthew Knepley 
1312000e7ae3SMatthew Knepley   Input Parameters:
1313000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1314000e7ae3SMatthew Knepley 
1315000e7ae3SMatthew Knepley   Level: developer
1316000e7ae3SMatthew Knepley 
1317000e7ae3SMatthew Knepley .keywords: TS, timestep
1318000e7ae3SMatthew Knepley @*/
1319000e7ae3SMatthew Knepley int TSDefaultPreStep(TS ts)
1320000e7ae3SMatthew Knepley {
1321000e7ae3SMatthew Knepley   PetscFunctionBegin;
1322000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1323000e7ae3SMatthew Knepley }
1324000e7ae3SMatthew Knepley 
1325e74ef692SMatthew Knepley #undef __FUNCT__
1326e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1327ac226902SBarry Smith /*@C
1328000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1329000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1330000e7ae3SMatthew Knepley   the time step.
1331000e7ae3SMatthew Knepley 
1332000e7ae3SMatthew Knepley   Collective on TS
1333000e7ae3SMatthew Knepley 
1334000e7ae3SMatthew Knepley   Input Parameters:
1335000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1336000e7ae3SMatthew Knepley - func - The function
1337000e7ae3SMatthew Knepley 
1338000e7ae3SMatthew Knepley   Calling sequence of func:
1339000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1340000e7ae3SMatthew Knepley 
1341000e7ae3SMatthew Knepley + t   - The current time
1342000e7ae3SMatthew Knepley - dt  - The current time step
1343000e7ae3SMatthew Knepley 
1344000e7ae3SMatthew Knepley   Level: intermediate
1345000e7ae3SMatthew Knepley 
1346000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1347000e7ae3SMatthew Knepley @*/
134889cef881SMatthew Knepley int TSSetUpdate(TS ts, int (*func)(TS, PetscReal, PetscReal *))
1349000e7ae3SMatthew Knepley {
1350000e7ae3SMatthew Knepley   PetscFunctionBegin;
13514482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1352000e7ae3SMatthew Knepley   ts->ops->update = func;
1353000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1354000e7ae3SMatthew Knepley }
1355000e7ae3SMatthew Knepley 
1356e74ef692SMatthew Knepley #undef __FUNCT__
1357e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1358000e7ae3SMatthew Knepley /*@
1359000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1360000e7ae3SMatthew Knepley 
1361000e7ae3SMatthew Knepley   Collective on TS
1362000e7ae3SMatthew Knepley 
1363000e7ae3SMatthew Knepley   Input Parameters:
1364000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1365000e7ae3SMatthew Knepley - t   - The current time
1366000e7ae3SMatthew Knepley 
1367000e7ae3SMatthew Knepley   Output Parameters:
1368000e7ae3SMatthew Knepley . dt  - The current time step
1369000e7ae3SMatthew Knepley 
1370000e7ae3SMatthew Knepley   Level: developer
1371000e7ae3SMatthew Knepley 
1372000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1373000e7ae3SMatthew Knepley @*/
137489cef881SMatthew Knepley int TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1375000e7ae3SMatthew Knepley {
1376000e7ae3SMatthew Knepley   PetscFunctionBegin;
1377000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1378000e7ae3SMatthew Knepley }
1379000e7ae3SMatthew Knepley 
1380e74ef692SMatthew Knepley #undef __FUNCT__
1381e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1382ac226902SBarry Smith /*@C
1383000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1384000e7ae3SMatthew Knepley   called once at the end of time stepping.
1385000e7ae3SMatthew Knepley 
1386000e7ae3SMatthew Knepley   Collective on TS
1387000e7ae3SMatthew Knepley 
1388000e7ae3SMatthew Knepley   Input Parameters:
1389000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1390000e7ae3SMatthew Knepley - func - The function
1391000e7ae3SMatthew Knepley 
1392000e7ae3SMatthew Knepley   Calling sequence of func:
1393000e7ae3SMatthew Knepley . func (TS ts);
1394000e7ae3SMatthew Knepley 
1395000e7ae3SMatthew Knepley   Level: intermediate
1396000e7ae3SMatthew Knepley 
1397000e7ae3SMatthew Knepley .keywords: TS, timestep
1398000e7ae3SMatthew Knepley @*/
1399000e7ae3SMatthew Knepley int TSSetPostStep(TS ts, int (*func)(TS))
1400000e7ae3SMatthew Knepley {
1401000e7ae3SMatthew Knepley   PetscFunctionBegin;
14024482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1403000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1404000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1405000e7ae3SMatthew Knepley }
1406000e7ae3SMatthew Knepley 
1407e74ef692SMatthew Knepley #undef __FUNCT__
1408e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1409000e7ae3SMatthew Knepley /*@
1410000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1411000e7ae3SMatthew Knepley 
1412000e7ae3SMatthew Knepley   Collective on TS
1413000e7ae3SMatthew Knepley 
1414000e7ae3SMatthew Knepley   Input Parameters:
1415000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1416000e7ae3SMatthew Knepley 
1417000e7ae3SMatthew Knepley   Level: developer
1418000e7ae3SMatthew Knepley 
1419000e7ae3SMatthew Knepley .keywords: TS, timestep
1420000e7ae3SMatthew Knepley @*/
1421000e7ae3SMatthew Knepley int TSDefaultPostStep(TS ts)
1422000e7ae3SMatthew Knepley {
1423000e7ae3SMatthew Knepley   PetscFunctionBegin;
1424000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1425000e7ae3SMatthew Knepley }
1426000e7ae3SMatthew Knepley 
1427d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1428d763cef2SBarry Smith 
14294a2ae208SSatish Balay #undef __FUNCT__
14304a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
1431d763cef2SBarry Smith /*@C
1432d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1433d763cef2SBarry Smith    timestep to display the iteration's  progress.
1434d763cef2SBarry Smith 
1435d763cef2SBarry Smith    Collective on TS
1436d763cef2SBarry Smith 
1437d763cef2SBarry Smith    Input Parameters:
1438d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1439d763cef2SBarry Smith .  func - monitoring routine
1440329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1441b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1442b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1443b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1444d763cef2SBarry Smith 
1445d763cef2SBarry Smith    Calling sequence of func:
144687828ca2SBarry Smith $    int func(TS ts,int steps,PetscReal time,Vec x,void *mctx)
1447d763cef2SBarry Smith 
1448d763cef2SBarry Smith +    ts - the TS context
1449d763cef2SBarry Smith .    steps - iteration number
14501f06c33eSBarry Smith .    time - current time
1451d763cef2SBarry Smith .    x - current iterate
1452d763cef2SBarry Smith -    mctx - [optional] monitoring context
1453d763cef2SBarry Smith 
1454d763cef2SBarry Smith    Notes:
1455d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1456d763cef2SBarry Smith    already has been loaded.
1457d763cef2SBarry Smith 
1458d763cef2SBarry Smith    Level: intermediate
1459d763cef2SBarry Smith 
1460d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1461d763cef2SBarry Smith 
1462d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
1463d763cef2SBarry Smith @*/
146487828ca2SBarry Smith int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*))
1465d763cef2SBarry Smith {
1466d763cef2SBarry Smith   PetscFunctionBegin;
14674482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1468d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
146929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1470d763cef2SBarry Smith   }
1471d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1472329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1473d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1474d763cef2SBarry Smith   PetscFunctionReturn(0);
1475d763cef2SBarry Smith }
1476d763cef2SBarry Smith 
14774a2ae208SSatish Balay #undef __FUNCT__
14784a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
1479d763cef2SBarry Smith /*@C
1480d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1481d763cef2SBarry Smith 
1482d763cef2SBarry Smith    Collective on TS
1483d763cef2SBarry Smith 
1484d763cef2SBarry Smith    Input Parameters:
1485d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1486d763cef2SBarry Smith 
1487d763cef2SBarry Smith    Notes:
1488d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1489d763cef2SBarry Smith 
1490d763cef2SBarry Smith    Level: intermediate
1491d763cef2SBarry Smith 
1492d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1493d763cef2SBarry Smith 
1494d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
1495d763cef2SBarry Smith @*/
1496d763cef2SBarry Smith int TSClearMonitor(TS ts)
1497d763cef2SBarry Smith {
1498d763cef2SBarry Smith   PetscFunctionBegin;
14994482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1500d763cef2SBarry Smith   ts->numbermonitors = 0;
1501d763cef2SBarry Smith   PetscFunctionReturn(0);
1502d763cef2SBarry Smith }
1503d763cef2SBarry Smith 
15044a2ae208SSatish Balay #undef __FUNCT__
15054a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
150687828ca2SBarry Smith int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
1507d763cef2SBarry Smith {
1508d132466eSBarry Smith   int ierr;
1509d132466eSBarry Smith 
1510d763cef2SBarry Smith   PetscFunctionBegin;
1511142b95e3SSatish Balay   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1512d763cef2SBarry Smith   PetscFunctionReturn(0);
1513d763cef2SBarry Smith }
1514d763cef2SBarry Smith 
15154a2ae208SSatish Balay #undef __FUNCT__
15164a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1517d763cef2SBarry Smith /*@
1518d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1519d763cef2SBarry Smith 
1520d763cef2SBarry Smith    Collective on TS
1521d763cef2SBarry Smith 
1522d763cef2SBarry Smith    Input Parameter:
1523d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1524d763cef2SBarry Smith 
1525d763cef2SBarry Smith    Output Parameters:
1526d763cef2SBarry Smith +  steps - number of iterations until termination
1527142b95e3SSatish Balay -  ptime - time until termination
1528d763cef2SBarry Smith 
1529d763cef2SBarry Smith    Level: beginner
1530d763cef2SBarry Smith 
1531d763cef2SBarry Smith .keywords: TS, timestep, solve
1532d763cef2SBarry Smith 
1533d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1534d763cef2SBarry Smith @*/
153587828ca2SBarry Smith int TSStep(TS ts,int *steps,PetscReal *ptime)
1536d763cef2SBarry Smith {
1537f1af5d2fSBarry Smith   int        ierr;
1538d763cef2SBarry Smith 
1539d763cef2SBarry Smith   PetscFunctionBegin;
15404482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1541d405a339SMatthew Knepley   if (!ts->setupcalled) {
1542d405a339SMatthew Knepley     ierr = TSSetUp(ts);                                                                                   CHKERRQ(ierr);
1543d405a339SMatthew Knepley   }
1544d405a339SMatthew Knepley 
1545d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);                                                        CHKERRQ(ierr);
1546d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);                                                                         CHKERRQ(ierr);
1547000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);                                                              CHKERRQ(ierr);
1548d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);                                                                        CHKERRQ(ierr);
1549d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);                                                          CHKERRQ(ierr);
1550d405a339SMatthew Knepley 
15514bb05414SBarry Smith   if (!PetscPreLoadingOn) {
15524bb05414SBarry Smith     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1553d405a339SMatthew Knepley   }
1554d763cef2SBarry Smith   PetscFunctionReturn(0);
1555d763cef2SBarry Smith }
1556d763cef2SBarry Smith 
15574a2ae208SSatish Balay #undef __FUNCT__
15584a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1559d763cef2SBarry Smith /*
1560d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1561d763cef2SBarry Smith */
156287828ca2SBarry Smith int TSMonitor(TS ts,int step,PetscReal ptime,Vec x)
1563d763cef2SBarry Smith {
1564d763cef2SBarry Smith   int i,ierr,n = ts->numbermonitors;
1565d763cef2SBarry Smith 
1566d763cef2SBarry Smith   PetscFunctionBegin;
1567d763cef2SBarry Smith   for (i=0; i<n; i++) {
1568142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1569d763cef2SBarry Smith   }
1570d763cef2SBarry Smith   PetscFunctionReturn(0);
1571d763cef2SBarry Smith }
1572d763cef2SBarry Smith 
1573d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1574d763cef2SBarry Smith 
15754a2ae208SSatish Balay #undef __FUNCT__
15764a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1577d763cef2SBarry Smith /*@C
1578d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1579d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1580d763cef2SBarry Smith 
1581d763cef2SBarry Smith    Collective on TS
1582d763cef2SBarry Smith 
1583d763cef2SBarry Smith    Input Parameters:
1584d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1585d763cef2SBarry Smith .  label - the title to put in the title bar
15867c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1587d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1588d763cef2SBarry Smith 
1589d763cef2SBarry Smith    Output Parameter:
1590d763cef2SBarry Smith .  draw - the drawing context
1591d763cef2SBarry Smith 
1592d763cef2SBarry Smith    Options Database Key:
1593d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1594d763cef2SBarry Smith 
1595d763cef2SBarry Smith    Notes:
1596b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1597d763cef2SBarry Smith 
1598d763cef2SBarry Smith    Level: intermediate
1599d763cef2SBarry Smith 
16007c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1601d763cef2SBarry Smith 
1602d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
16037c922b88SBarry Smith 
1604d763cef2SBarry Smith @*/
16051836bdbcSSatish Balay int TSLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1606d763cef2SBarry Smith {
1607b0a32e0cSBarry Smith   PetscDraw win;
1608d763cef2SBarry Smith   int       ierr;
1609d763cef2SBarry Smith 
1610d763cef2SBarry Smith   PetscFunctionBegin;
1611b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1612b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1613b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1614b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1615d763cef2SBarry Smith 
1616b0a32e0cSBarry Smith   PetscLogObjectParent(*draw,win);
1617d763cef2SBarry Smith   PetscFunctionReturn(0);
1618d763cef2SBarry Smith }
1619d763cef2SBarry Smith 
16204a2ae208SSatish Balay #undef __FUNCT__
16214a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
162287828ca2SBarry Smith int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx)
1623d763cef2SBarry Smith {
1624b0a32e0cSBarry Smith   PetscDrawLG lg = (PetscDrawLG) monctx;
162587828ca2SBarry Smith   PetscReal      x,y = ptime;
1626d763cef2SBarry Smith   int         ierr;
1627d763cef2SBarry Smith 
1628d763cef2SBarry Smith   PetscFunctionBegin;
16297c922b88SBarry Smith   if (!monctx) {
16307c922b88SBarry Smith     MPI_Comm    comm;
1631b0a32e0cSBarry Smith     PetscViewer viewer;
16327c922b88SBarry Smith 
16337c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1634b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1635b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
16367c922b88SBarry Smith   }
16377c922b88SBarry Smith 
1638b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
163987828ca2SBarry Smith   x = (PetscReal)n;
1640b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1641d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1642b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1643d763cef2SBarry Smith   }
1644d763cef2SBarry Smith   PetscFunctionReturn(0);
1645d763cef2SBarry Smith }
1646d763cef2SBarry Smith 
16474a2ae208SSatish Balay #undef __FUNCT__
16484a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1649d763cef2SBarry Smith /*@C
1650d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1651d763cef2SBarry Smith    with TSLGMonitorCreate().
1652d763cef2SBarry Smith 
1653b0a32e0cSBarry Smith    Collective on PetscDrawLG
1654d763cef2SBarry Smith 
1655d763cef2SBarry Smith    Input Parameter:
1656d763cef2SBarry Smith .  draw - the drawing context
1657d763cef2SBarry Smith 
1658d763cef2SBarry Smith    Level: intermediate
1659d763cef2SBarry Smith 
1660d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1661d763cef2SBarry Smith 
1662d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1663d763cef2SBarry Smith @*/
1664b0a32e0cSBarry Smith int TSLGMonitorDestroy(PetscDrawLG drawlg)
1665d763cef2SBarry Smith {
1666b0a32e0cSBarry Smith   PetscDraw draw;
1667d763cef2SBarry Smith   int       ierr;
1668d763cef2SBarry Smith 
1669d763cef2SBarry Smith   PetscFunctionBegin;
1670b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1671b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1672b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1673d763cef2SBarry Smith   PetscFunctionReturn(0);
1674d763cef2SBarry Smith }
1675d763cef2SBarry Smith 
16764a2ae208SSatish Balay #undef __FUNCT__
16774a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1678d763cef2SBarry Smith /*@
1679d763cef2SBarry Smith    TSGetTime - Gets the current time.
1680d763cef2SBarry Smith 
1681d763cef2SBarry Smith    Not Collective
1682d763cef2SBarry Smith 
1683d763cef2SBarry Smith    Input Parameter:
1684d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1685d763cef2SBarry Smith 
1686d763cef2SBarry Smith    Output Parameter:
1687d763cef2SBarry Smith .  t  - the current time
1688d763cef2SBarry Smith 
1689d763cef2SBarry Smith    Contributed by: Matthew Knepley
1690d763cef2SBarry Smith 
1691d763cef2SBarry Smith    Level: beginner
1692d763cef2SBarry Smith 
1693d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1694d763cef2SBarry Smith 
1695d763cef2SBarry Smith .keywords: TS, get, time
1696d763cef2SBarry Smith @*/
169787828ca2SBarry Smith int TSGetTime(TS ts,PetscReal* t)
1698d763cef2SBarry Smith {
1699d763cef2SBarry Smith   PetscFunctionBegin;
17004482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
17014482741eSBarry Smith   PetscValidDoublePointer(t,2);
1702d763cef2SBarry Smith   *t = ts->ptime;
1703d763cef2SBarry Smith   PetscFunctionReturn(0);
1704d763cef2SBarry Smith }
1705d763cef2SBarry Smith 
17064a2ae208SSatish Balay #undef __FUNCT__
17074a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1708d763cef2SBarry Smith /*@C
1709d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1710d763cef2SBarry Smith    TS options in the database.
1711d763cef2SBarry Smith 
1712d763cef2SBarry Smith    Collective on TS
1713d763cef2SBarry Smith 
1714d763cef2SBarry Smith    Input Parameter:
1715d763cef2SBarry Smith +  ts     - The TS context
1716d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1717d763cef2SBarry Smith 
1718d763cef2SBarry Smith    Notes:
1719d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1720d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1721d763cef2SBarry Smith    hyphen.
1722d763cef2SBarry Smith 
1723d763cef2SBarry Smith    Contributed by: Matthew Knepley
1724d763cef2SBarry Smith 
1725d763cef2SBarry Smith    Level: advanced
1726d763cef2SBarry Smith 
1727d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1728d763cef2SBarry Smith 
1729d763cef2SBarry Smith .seealso: TSSetFromOptions()
1730d763cef2SBarry Smith 
1731d763cef2SBarry Smith @*/
17321836bdbcSSatish Balay int TSSetOptionsPrefix(TS ts,const char prefix[])
1733d763cef2SBarry Smith {
1734d763cef2SBarry Smith   int ierr;
1735d763cef2SBarry Smith 
1736d763cef2SBarry Smith   PetscFunctionBegin;
17374482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1738d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1739d763cef2SBarry Smith   switch(ts->problem_type) {
1740d763cef2SBarry Smith     case TS_NONLINEAR:
1741d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1742d763cef2SBarry Smith       break;
1743d763cef2SBarry Smith     case TS_LINEAR:
174494b7f48cSBarry Smith       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1745d763cef2SBarry Smith       break;
1746d763cef2SBarry Smith   }
1747d763cef2SBarry Smith   PetscFunctionReturn(0);
1748d763cef2SBarry Smith }
1749d763cef2SBarry Smith 
1750d763cef2SBarry Smith 
17514a2ae208SSatish Balay #undef __FUNCT__
17524a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1753d763cef2SBarry Smith /*@C
1754d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1755d763cef2SBarry Smith    TS options in the database.
1756d763cef2SBarry Smith 
1757d763cef2SBarry Smith    Collective on TS
1758d763cef2SBarry Smith 
1759d763cef2SBarry Smith    Input Parameter:
1760d763cef2SBarry Smith +  ts     - The TS context
1761d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1762d763cef2SBarry Smith 
1763d763cef2SBarry Smith    Notes:
1764d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1765d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1766d763cef2SBarry Smith    hyphen.
1767d763cef2SBarry Smith 
1768d763cef2SBarry Smith    Contributed by: Matthew Knepley
1769d763cef2SBarry Smith 
1770d763cef2SBarry Smith    Level: advanced
1771d763cef2SBarry Smith 
1772d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1773d763cef2SBarry Smith 
1774d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1775d763cef2SBarry Smith 
1776d763cef2SBarry Smith @*/
17771836bdbcSSatish Balay int TSAppendOptionsPrefix(TS ts,const char prefix[])
1778d763cef2SBarry Smith {
1779d763cef2SBarry Smith   int ierr;
1780d763cef2SBarry Smith 
1781d763cef2SBarry Smith   PetscFunctionBegin;
17824482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1783d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1784d763cef2SBarry Smith   switch(ts->problem_type) {
1785d763cef2SBarry Smith     case TS_NONLINEAR:
1786d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1787d763cef2SBarry Smith       break;
1788d763cef2SBarry Smith     case TS_LINEAR:
178994b7f48cSBarry Smith       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1790d763cef2SBarry Smith       break;
1791d763cef2SBarry Smith   }
1792d763cef2SBarry Smith   PetscFunctionReturn(0);
1793d763cef2SBarry Smith }
1794d763cef2SBarry Smith 
17954a2ae208SSatish Balay #undef __FUNCT__
17964a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1797d763cef2SBarry Smith /*@C
1798d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1799d763cef2SBarry Smith    TS options in the database.
1800d763cef2SBarry Smith 
1801d763cef2SBarry Smith    Not Collective
1802d763cef2SBarry Smith 
1803d763cef2SBarry Smith    Input Parameter:
1804d763cef2SBarry Smith .  ts - The TS context
1805d763cef2SBarry Smith 
1806d763cef2SBarry Smith    Output Parameter:
1807d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1808d763cef2SBarry Smith 
1809d763cef2SBarry Smith    Contributed by: Matthew Knepley
1810d763cef2SBarry Smith 
1811d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1812d763cef2SBarry Smith    sufficient length to hold the prefix.
1813d763cef2SBarry Smith 
1814d763cef2SBarry Smith    Level: intermediate
1815d763cef2SBarry Smith 
1816d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1817d763cef2SBarry Smith 
1818d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1819d763cef2SBarry Smith @*/
18201836bdbcSSatish Balay int TSGetOptionsPrefix(TS ts,char *prefix[])
1821d763cef2SBarry Smith {
1822d763cef2SBarry Smith   int ierr;
1823d763cef2SBarry Smith 
1824d763cef2SBarry Smith   PetscFunctionBegin;
18254482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
18264482741eSBarry Smith   PetscValidPointer(prefix,2);
1827d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1828d763cef2SBarry Smith   PetscFunctionReturn(0);
1829d763cef2SBarry Smith }
1830d763cef2SBarry Smith 
18314a2ae208SSatish Balay #undef __FUNCT__
18324a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1833d763cef2SBarry Smith /*@C
1834d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1835d763cef2SBarry Smith 
1836d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1837d763cef2SBarry Smith 
1838d763cef2SBarry Smith    Input Parameter:
1839d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1840d763cef2SBarry Smith 
1841d763cef2SBarry Smith    Output Parameters:
1842d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1843d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1844d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1845d763cef2SBarry Smith 
1846d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1847d763cef2SBarry Smith 
1848d763cef2SBarry Smith    Contributed by: Matthew Knepley
1849d763cef2SBarry Smith 
1850d763cef2SBarry Smith    Level: intermediate
1851d763cef2SBarry Smith 
1852d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1853d763cef2SBarry Smith 
1854d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1855d763cef2SBarry Smith 
1856d763cef2SBarry Smith @*/
1857d763cef2SBarry Smith int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1858d763cef2SBarry Smith {
1859d763cef2SBarry Smith   PetscFunctionBegin;
18604482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1861d763cef2SBarry Smith   if (A)   *A = ts->A;
1862d763cef2SBarry Smith   if (M)   *M = ts->B;
1863d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1864d763cef2SBarry Smith   PetscFunctionReturn(0);
1865d763cef2SBarry Smith }
1866d763cef2SBarry Smith 
18674a2ae208SSatish Balay #undef __FUNCT__
18684a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1869d763cef2SBarry Smith /*@C
1870d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1871d763cef2SBarry Smith 
1872d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1873d763cef2SBarry Smith 
1874d763cef2SBarry Smith    Input Parameter:
1875d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1876d763cef2SBarry Smith 
1877d763cef2SBarry Smith    Output Parameters:
1878d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1879d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1880d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1881d763cef2SBarry Smith 
1882d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1883d763cef2SBarry Smith 
1884d763cef2SBarry Smith    Contributed by: Matthew Knepley
1885d763cef2SBarry Smith 
1886d763cef2SBarry Smith    Level: intermediate
1887d763cef2SBarry Smith 
1888d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1889d763cef2SBarry Smith 
1890d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1891d763cef2SBarry Smith @*/
1892d763cef2SBarry Smith int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1893d763cef2SBarry Smith {
1894d763cef2SBarry Smith   int ierr;
1895d763cef2SBarry Smith 
1896d763cef2SBarry Smith   PetscFunctionBegin;
1897d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1898d763cef2SBarry Smith   PetscFunctionReturn(0);
1899d763cef2SBarry Smith }
1900d763cef2SBarry Smith 
19011713a123SBarry Smith #undef __FUNCT__
19021713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
19031713a123SBarry Smith /*@C
19041713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
19051713a123SBarry Smith    VecView() for the solution at each timestep
19061713a123SBarry Smith 
19071713a123SBarry Smith    Collective on TS
19081713a123SBarry Smith 
19091713a123SBarry Smith    Input Parameters:
19101713a123SBarry Smith +  ts - the TS context
19111713a123SBarry Smith .  step - current time-step
1912142b95e3SSatish Balay .  ptime - current time
19131713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
19141713a123SBarry Smith 
19151713a123SBarry Smith    Level: intermediate
19161713a123SBarry Smith 
19171713a123SBarry Smith .keywords: TS,  vector, monitor, view
19181713a123SBarry Smith 
19191713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
19201713a123SBarry Smith @*/
1921142b95e3SSatish Balay int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
19221713a123SBarry Smith {
19231713a123SBarry Smith   int         ierr;
19241713a123SBarry Smith   PetscViewer viewer = (PetscViewer) dummy;
19251713a123SBarry Smith 
19261713a123SBarry Smith   PetscFunctionBegin;
19271713a123SBarry Smith   if (!viewer) {
19281713a123SBarry Smith     MPI_Comm comm;
19291713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
19301713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
19311713a123SBarry Smith   }
19321713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
19331713a123SBarry Smith   PetscFunctionReturn(0);
19341713a123SBarry Smith }
19351713a123SBarry Smith 
19361713a123SBarry Smith 
19371713a123SBarry Smith 
1938