xref: /petsc/src/ts/interface/ts.c (revision 95f0b56261087137261fecb6c5c0f94357e5c2d7)
163dd3a1aSKris Buschelman #define PETSCTS_DLL
263dd3a1aSKris Buschelman 
3b9147fbbSdalcinl #include "include/private/tsimpl.h"        /*I "petscts.h"  I*/
4d763cef2SBarry Smith 
5d5ba7fb7SMatthew Knepley /* Logging support */
663dd3a1aSKris Buschelman PetscCookie PETSCTS_DLLEXPORT TS_COOKIE = 0;
76849ba73SBarry Smith PetscEvent  TS_Step = 0, TS_PseudoComputeTimeStep = 0, TS_FunctionEval = 0, TS_JacobianEval = 0;
8d405a339SMatthew Knepley 
94a2ae208SSatish Balay #undef __FUNCT__
10bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
11bdad233fSMatthew Knepley /*
12bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
13bdad233fSMatthew Knepley 
14bdad233fSMatthew Knepley   Collective on TS
15bdad233fSMatthew Knepley 
16bdad233fSMatthew Knepley   Input Parameter:
17bdad233fSMatthew Knepley . ts - The ts
18bdad233fSMatthew Knepley 
19bdad233fSMatthew Knepley   Level: intermediate
20bdad233fSMatthew Knepley 
21bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
22bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
23bdad233fSMatthew Knepley */
246849ba73SBarry Smith static PetscErrorCode TSSetTypeFromOptions(TS ts)
25bdad233fSMatthew Knepley {
26bdad233fSMatthew Knepley   PetscTruth     opt;
272fc52814SBarry Smith   const char     *defaultType;
28bdad233fSMatthew Knepley   char           typeName[256];
29dfbe8321SBarry Smith   PetscErrorCode ierr;
30bdad233fSMatthew Knepley 
31bdad233fSMatthew Knepley   PetscFunctionBegin;
32abc0a331SBarry Smith   if (ts->type_name) {
33bdad233fSMatthew Knepley     defaultType = ts->type_name;
34bdad233fSMatthew Knepley   } else {
35bdad233fSMatthew Knepley     defaultType = TS_EULER;
36bdad233fSMatthew Knepley   }
37bdad233fSMatthew Knepley 
38bdad233fSMatthew Knepley   if (!TSRegisterAllCalled) {
39bdad233fSMatthew Knepley     ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);
40bdad233fSMatthew Knepley   }
41bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
42a7cc72afSBarry Smith   if (opt) {
43bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
44bdad233fSMatthew Knepley   } else {
45bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
46bdad233fSMatthew Knepley   }
47bdad233fSMatthew Knepley   PetscFunctionReturn(0);
48bdad233fSMatthew Knepley }
49bdad233fSMatthew Knepley 
50bdad233fSMatthew Knepley #undef __FUNCT__
51bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
52bdad233fSMatthew Knepley /*@
53bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
54bdad233fSMatthew Knepley 
55bdad233fSMatthew Knepley    Collective on TS
56bdad233fSMatthew Knepley 
57bdad233fSMatthew Knepley    Input Parameter:
58bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
59bdad233fSMatthew Knepley 
60bdad233fSMatthew Knepley    Options Database Keys:
610f3b3ca1SHong Zhang +  -ts_type <type> - TS_EULER, TS_BEULER, TS_SUNDIALS, TS_PSEUDO, TS_CRANK_NICHOLSON
62bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
63bdad233fSMatthew Knepley .  -ts_max_time time - maximum time to compute to
64bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
65bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
66a6570f20SBarry Smith -  -ts_monitor_draw - plot information at each timestep
67bdad233fSMatthew Knepley 
68bdad233fSMatthew Knepley    Level: beginner
69bdad233fSMatthew Knepley 
70bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
71bdad233fSMatthew Knepley 
72bdad233fSMatthew Knepley .seealso: TSGetType
73bdad233fSMatthew Knepley @*/
7463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetFromOptions(TS ts)
75bdad233fSMatthew Knepley {
76bdad233fSMatthew Knepley   PetscReal               dt;
77eabae89aSBarry Smith   PetscTruth              opt,flg;
78dfbe8321SBarry Smith   PetscErrorCode          ierr;
79a34d58ebSBarry Smith   PetscViewerASCIIMonitor monviewer;
80eabae89aSBarry Smith   char                    monfilename[PETSC_MAX_PATH_LEN];
81bdad233fSMatthew Knepley 
82bdad233fSMatthew Knepley   PetscFunctionBegin;
834482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
84bdad233fSMatthew Knepley   ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");CHKERRQ(ierr);
85bdad233fSMatthew Knepley 
86bdad233fSMatthew Knepley     /* Handle generic TS options */
87bdad233fSMatthew Knepley     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
88bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
89bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
90bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
91a7cc72afSBarry Smith     if (opt) {
92bdad233fSMatthew Knepley       ts->initial_time_step = ts->time_step = dt;
93bdad233fSMatthew Knepley     }
94bdad233fSMatthew Knepley 
95bdad233fSMatthew Knepley     /* Monitor options */
96a6570f20SBarry Smith     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
97eabae89aSBarry Smith     if (flg) {
98a34d58ebSBarry Smith       ierr = PetscViewerASCIIMonitorCreate(ts->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
99a34d58ebSBarry Smith       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
100bdad233fSMatthew Knepley     }
101a6570f20SBarry Smith     ierr = PetscOptionsName("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",&opt);CHKERRQ(ierr);
102a7cc72afSBarry Smith     if (opt) {
103a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
104bdad233fSMatthew Knepley     }
105a6570f20SBarry Smith     ierr = PetscOptionsName("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",&opt);CHKERRQ(ierr);
106a7cc72afSBarry Smith     if (opt) {
107a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
108bdad233fSMatthew Knepley     }
109bdad233fSMatthew Knepley 
110bdad233fSMatthew Knepley     /* Handle TS type options */
111bdad233fSMatthew Knepley     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
112bdad233fSMatthew Knepley 
113bdad233fSMatthew Knepley     /* Handle specific TS options */
114abc0a331SBarry Smith     if (ts->ops->setfromoptions) {
115bdad233fSMatthew Knepley       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
116bdad233fSMatthew Knepley     }
117bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
118bdad233fSMatthew Knepley 
119bdad233fSMatthew Knepley   /* Handle subobject options */
120bdad233fSMatthew Knepley   switch(ts->problem_type) {
121156fc9a6SMatthew Knepley     /* Should check for implicit/explicit */
122bdad233fSMatthew Knepley   case TS_LINEAR:
123abc0a331SBarry Smith     if (ts->ksp) {
1248beb423aSHong Zhang       ierr = KSPSetOperators(ts->ksp,ts->Arhs,ts->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
12594b7f48cSBarry Smith       ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr);
126156fc9a6SMatthew Knepley     }
127bdad233fSMatthew Knepley     break;
128bdad233fSMatthew Knepley   case TS_NONLINEAR:
129abc0a331SBarry Smith     if (ts->snes) {
1307c236d22SBarry Smith       /* this is a bit of a hack, but it gets the matrix information into SNES earlier
1317c236d22SBarry Smith          so that SNES and KSP have more information to pick reasonable defaults
1327c236d22SBarry Smith          before they allow users to set options */
1338beb423aSHong Zhang       ierr = SNESSetJacobian(ts->snes,ts->Arhs,ts->B,0,ts);CHKERRQ(ierr);
134bdad233fSMatthew Knepley       ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
135156fc9a6SMatthew Knepley     }
136bdad233fSMatthew Knepley     break;
137bdad233fSMatthew Knepley   default:
13877431f27SBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type);
139bdad233fSMatthew Knepley   }
140bdad233fSMatthew Knepley 
141bdad233fSMatthew Knepley   PetscFunctionReturn(0);
142bdad233fSMatthew Knepley }
143bdad233fSMatthew Knepley 
144bdad233fSMatthew Knepley #undef  __FUNCT__
145bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
146bdad233fSMatthew Knepley /*@
147bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
148bdad233fSMatthew Knepley 
149bdad233fSMatthew Knepley   Collective on TS
150bdad233fSMatthew Knepley 
151bdad233fSMatthew Knepley   Input Parameter:
152bdad233fSMatthew Knepley . ts - The ts
153bdad233fSMatthew Knepley 
154bdad233fSMatthew Knepley   Level: intermediate
155bdad233fSMatthew Knepley 
156bdad233fSMatthew Knepley .keywords: TS, view, options, database
157bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
158bdad233fSMatthew Knepley @*/
15963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSViewFromOptions(TS ts,const char title[])
160bdad233fSMatthew Knepley {
161bdad233fSMatthew Knepley   PetscViewer    viewer;
162bdad233fSMatthew Knepley   PetscDraw      draw;
163bdad233fSMatthew Knepley   PetscTruth     opt;
164e10c49a3SBarry Smith   char           fileName[PETSC_MAX_PATH_LEN];
165dfbe8321SBarry Smith   PetscErrorCode ierr;
166bdad233fSMatthew Knepley 
167bdad233fSMatthew Knepley   PetscFunctionBegin;
168eabae89aSBarry Smith   ierr = PetscOptionsGetString(ts->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr);
169eabae89aSBarry Smith   if (opt && !PetscPreLoadingOn) {
170eabae89aSBarry Smith     ierr = PetscViewerASCIIOpen(ts->comm,fileName,&viewer);CHKERRQ(ierr);
171bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
172bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
173bdad233fSMatthew Knepley   }
174bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);CHKERRQ(ierr);
175a7cc72afSBarry Smith   if (opt) {
176bdad233fSMatthew Knepley     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
177bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
178a7cc72afSBarry Smith     if (title) {
1791836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
180bdad233fSMatthew Knepley     } else {
181bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject) ts);CHKERRQ(ierr);
1821836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, ts->name);CHKERRQ(ierr);
183bdad233fSMatthew Knepley     }
184bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
185bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
186bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
187bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
188bdad233fSMatthew Knepley   }
189bdad233fSMatthew Knepley   PetscFunctionReturn(0);
190bdad233fSMatthew Knepley }
191bdad233fSMatthew Knepley 
192bdad233fSMatthew Knepley #undef __FUNCT__
1934a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
194a7a1495cSBarry Smith /*@
1958c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
196a7a1495cSBarry Smith       set with TSSetRHSJacobian().
197a7a1495cSBarry Smith 
198a7a1495cSBarry Smith    Collective on TS and Vec
199a7a1495cSBarry Smith 
200a7a1495cSBarry Smith    Input Parameters:
201a7a1495cSBarry Smith +  ts - the SNES context
202a7a1495cSBarry Smith .  t - current timestep
203a7a1495cSBarry Smith -  x - input vector
204a7a1495cSBarry Smith 
205a7a1495cSBarry Smith    Output Parameters:
206a7a1495cSBarry Smith +  A - Jacobian matrix
207a7a1495cSBarry Smith .  B - optional preconditioning matrix
208a7a1495cSBarry Smith -  flag - flag indicating matrix structure
209a7a1495cSBarry Smith 
210a7a1495cSBarry Smith    Notes:
211a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
212a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
213a7a1495cSBarry Smith 
21494b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
215a7a1495cSBarry Smith    flag parameter.
216a7a1495cSBarry Smith 
217a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
218a7a1495cSBarry Smith 
219a7a1495cSBarry Smith    Level: developer
220a7a1495cSBarry Smith 
221a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
222a7a1495cSBarry Smith 
22394b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
224a7a1495cSBarry Smith @*/
22563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
226a7a1495cSBarry Smith {
227dfbe8321SBarry Smith   PetscErrorCode ierr;
228a7a1495cSBarry Smith 
229a7a1495cSBarry Smith   PetscFunctionBegin;
2304482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2314482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
232c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
233a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
23429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
235a7a1495cSBarry Smith   }
236000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
237d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
238a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
239a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
240000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
241a7a1495cSBarry Smith     PetscStackPop;
242d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
243a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2444482741eSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE,4);
2454482741eSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE,5);
246ef66eb69SBarry Smith   } else {
247ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249ef66eb69SBarry Smith     if (*A != *B) {
250ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252ef66eb69SBarry Smith     }
253ef66eb69SBarry Smith   }
254a7a1495cSBarry Smith   PetscFunctionReturn(0);
255a7a1495cSBarry Smith }
256a7a1495cSBarry Smith 
2574a2ae208SSatish Balay #undef __FUNCT__
2584a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
259d763cef2SBarry Smith /*
260d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
261d763cef2SBarry Smith 
262d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
263d763cef2SBarry Smith    this routine applies the matrix.
264d763cef2SBarry Smith */
265dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
266d763cef2SBarry Smith {
267dfbe8321SBarry Smith   PetscErrorCode ierr;
268d763cef2SBarry Smith 
269d763cef2SBarry Smith   PetscFunctionBegin;
2704482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2714482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
2724482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
273d763cef2SBarry Smith 
274d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
275000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
276d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
277000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
278d763cef2SBarry Smith     PetscStackPop;
2797c922b88SBarry Smith   } else {
280000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
281d763cef2SBarry Smith       MatStructure flg;
282d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
2838beb423aSHong Zhang       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
284d763cef2SBarry Smith       PetscStackPop;
285d763cef2SBarry Smith     }
2868beb423aSHong Zhang     ierr = MatMult(ts->Arhs,x,y);CHKERRQ(ierr);
2877c922b88SBarry Smith   }
288d763cef2SBarry Smith 
289d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
290d763cef2SBarry Smith 
291d763cef2SBarry Smith   PetscFunctionReturn(0);
292d763cef2SBarry Smith }
293d763cef2SBarry Smith 
2944a2ae208SSatish Balay #undef __FUNCT__
2954a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
296d763cef2SBarry Smith /*@C
297d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
298d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
299d763cef2SBarry Smith 
300d763cef2SBarry Smith     Collective on TS
301d763cef2SBarry Smith 
302d763cef2SBarry Smith     Input Parameters:
303d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
304d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
305d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
306d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
307d763cef2SBarry Smith 
308d763cef2SBarry Smith     Calling sequence of func:
30987828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
310d763cef2SBarry Smith 
311d763cef2SBarry Smith +   t - current timestep
312d763cef2SBarry Smith .   u - input vector
313d763cef2SBarry Smith .   F - function vector
314d763cef2SBarry Smith -   ctx - [optional] user-defined function context
315d763cef2SBarry Smith 
316d763cef2SBarry Smith     Important:
317d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
318d763cef2SBarry Smith 
319d763cef2SBarry Smith     Level: beginner
320d763cef2SBarry Smith 
321d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
322d763cef2SBarry Smith 
323d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
324d763cef2SBarry Smith @*/
32563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
326d763cef2SBarry Smith {
327d763cef2SBarry Smith   PetscFunctionBegin;
328d763cef2SBarry Smith 
3294482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
330d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
33129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
332d763cef2SBarry Smith   }
333000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
334d763cef2SBarry Smith   ts->funP             = ctx;
335d763cef2SBarry Smith   PetscFunctionReturn(0);
336d763cef2SBarry Smith }
337d763cef2SBarry Smith 
3384a2ae208SSatish Balay #undef __FUNCT__
339*95f0b562SHong Zhang #define __FUNCT__ "TSSetMatrices"
340*95f0b562SHong Zhang /*@C
341*95f0b562SHong Zhang    TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs,
342*95f0b562SHong Zhang    where Alhs(t) U_t = Arhs(t) U.
343*95f0b562SHong Zhang 
344*95f0b562SHong Zhang    Collective on TS
345*95f0b562SHong Zhang 
346*95f0b562SHong Zhang    Input Parameters:
347*95f0b562SHong Zhang +  ts   - the TS context obtained from TSCreate()
348*95f0b562SHong Zhang .  Arhs - matrix
349*95f0b562SHong Zhang .  frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
350*95f0b562SHong Zhang           if Arhs is not a function of t.
351*95f0b562SHong Zhang .  Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix.
352*95f0b562SHong Zhang .  flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
353*95f0b562SHong Zhang           if Alhs is not a function of t.
354*95f0b562SHong Zhang .  flag - flag indicating information about the matrix structure of Arhs and Alhs.
355*95f0b562SHong Zhang           The available options are
356*95f0b562SHong Zhang             SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs
357*95f0b562SHong Zhang             DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs
358*95f0b562SHong Zhang -  ctx  - [optional] user-defined context for private data for the
359*95f0b562SHong Zhang           matrix evaluation routine (may be PETSC_NULL)
360*95f0b562SHong Zhang 
361*95f0b562SHong Zhang    Calling sequence of func:
362*95f0b562SHong Zhang $     func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
363*95f0b562SHong Zhang 
364*95f0b562SHong Zhang +  t - current timestep
365*95f0b562SHong Zhang .  A - matrix A, where U_t = A(t) U
366*95f0b562SHong Zhang .  B - preconditioner matrix, usually the same as A
367*95f0b562SHong Zhang .  flag - flag indicating information about the preconditioner matrix
368*95f0b562SHong Zhang           structure (same as flag in KSPSetOperators())
369*95f0b562SHong Zhang -  ctx - [optional] user-defined context for matrix evaluation routine
370*95f0b562SHong Zhang 
371*95f0b562SHong Zhang    Notes:
372*95f0b562SHong Zhang    The routine func() takes Mat* as the matrix arguments rather than Mat.
373*95f0b562SHong Zhang    This allows the matrix evaluation routine to replace Arhs or Alhs with a
374*95f0b562SHong Zhang    completely new new matrix structure (not just different matrix elements)
375*95f0b562SHong Zhang    when appropriate, for instance, if the nonzero structure is changing
376*95f0b562SHong Zhang    throughout the global iterations.
377*95f0b562SHong Zhang 
378*95f0b562SHong Zhang    Important:
379*95f0b562SHong Zhang    The user MUST call either this routine or TSSetRHSFunction().
380*95f0b562SHong Zhang 
381*95f0b562SHong Zhang    Level: beginner
382*95f0b562SHong Zhang 
383*95f0b562SHong Zhang .keywords: TS, timestep, set, matrix
384*95f0b562SHong Zhang 
385*95f0b562SHong Zhang .seealso: TSSetRHSFunction()
386*95f0b562SHong Zhang @*/
387*95f0b562SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSetMatrices(TS ts,Mat Arhs,PetscErrorCode (*frhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat Alhs,PetscErrorCode (*flhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure flag,void *ctx)
388*95f0b562SHong Zhang {
389*95f0b562SHong Zhang   PetscFunctionBegin;
390*95f0b562SHong Zhang   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
391*95f0b562SHong Zhang   PetscValidHeaderSpecific(Arhs,MAT_COOKIE,2);
392*95f0b562SHong Zhang   PetscCheckSameComm(ts,1,Arhs,2);
393*95f0b562SHong Zhang   if (Alhs) PetscCheckSameComm(ts,1,Alhs,4);
394*95f0b562SHong Zhang   if (ts->problem_type == TS_NONLINEAR)
395*95f0b562SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
396*95f0b562SHong Zhang 
397*95f0b562SHong Zhang   ts->ops->rhsmatrix = frhs;
398*95f0b562SHong Zhang   ts->ops->lhsmatrix = flhs;
399*95f0b562SHong Zhang   ts->jacP           = ctx;
400*95f0b562SHong Zhang   ts->Arhs           = Arhs;
401*95f0b562SHong Zhang   ts->Alhs           = Alhs;
402*95f0b562SHong Zhang   ts->matflg         = flag;
403*95f0b562SHong Zhang   PetscFunctionReturn(0);
404*95f0b562SHong Zhang }
405*95f0b562SHong Zhang 
406*95f0b562SHong Zhang #undef __FUNCT__
4074a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
408d763cef2SBarry Smith /*@C
409d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
410d763cef2SBarry Smith    Also sets the location to store A.
411d763cef2SBarry Smith 
412d763cef2SBarry Smith    Collective on TS
413d763cef2SBarry Smith 
414d763cef2SBarry Smith    Input Parameters:
415d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
416d763cef2SBarry Smith .  A   - matrix
417d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
418d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
419d763cef2SBarry Smith          if A is not a function of t.
420d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
421d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
422d763cef2SBarry Smith 
423d763cef2SBarry Smith    Calling sequence of func:
424a7cc72afSBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
425d763cef2SBarry Smith 
426d763cef2SBarry Smith +  t - current timestep
427d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
428d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
429d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
43094b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
431d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
432d763cef2SBarry Smith 
433d763cef2SBarry Smith    Notes:
43494b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
435d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
436d763cef2SBarry Smith 
437d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
438d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
439d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
440d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
441d763cef2SBarry Smith    throughout the global iterations.
442d763cef2SBarry Smith 
443d763cef2SBarry Smith    Important:
444d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
445d763cef2SBarry Smith 
446d763cef2SBarry Smith    Level: beginner
447d763cef2SBarry Smith 
448d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
449d763cef2SBarry Smith 
450d763cef2SBarry Smith .seealso: TSSetRHSFunction()
451d763cef2SBarry Smith @*/
45263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSMatrix(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
453d763cef2SBarry Smith {
454d763cef2SBarry Smith   PetscFunctionBegin;
4554482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
4564482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
4574482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
458c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
459c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
460d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
46129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
462d763cef2SBarry Smith   }
463d763cef2SBarry Smith 
464000e7ae3SMatthew Knepley   ts->ops->rhsmatrix = f;
465d763cef2SBarry Smith   ts->jacP           = ctx;
4668beb423aSHong Zhang   ts->Arhs           = A;
467d763cef2SBarry Smith   ts->B              = B;
468d763cef2SBarry Smith   PetscFunctionReturn(0);
469d763cef2SBarry Smith }
470d763cef2SBarry Smith 
4714a2ae208SSatish Balay #undef __FUNCT__
472aa644b49SHong Zhang #define __FUNCT__ "TSSetLHSMatrix"
473aa644b49SHong Zhang /*@C
474aa644b49SHong Zhang    TSSetLHSMatrix - Sets the function to compute the matrix A, where A U_t = F(U).
475aa644b49SHong Zhang    Also sets the location to store A.
476aa644b49SHong Zhang 
477aa644b49SHong Zhang    Collective on TS
478aa644b49SHong Zhang 
479aa644b49SHong Zhang    Input Parameters:
480aa644b49SHong Zhang +  ts  - the TS context obtained from TSCreate()
481aa644b49SHong Zhang .  A   - matrix
4828beb423aSHong Zhang .  B   - ignored
483aa644b49SHong Zhang .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
484aa644b49SHong Zhang          if A is not a function of t.
485aa644b49SHong Zhang -  ctx - [optional] user-defined context for private data for the
486aa644b49SHong Zhang           matrix evaluation routine (may be PETSC_NULL)
487aa644b49SHong Zhang 
488aa644b49SHong Zhang    Calling sequence of func:
489aa644b49SHong Zhang $     func (TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
490aa644b49SHong Zhang 
491aa644b49SHong Zhang +  t - current timestep
492aa644b49SHong Zhang .  A - matrix A, where A U_t = F(U)
4938beb423aSHong Zhang .  B - ignored
494aa644b49SHong Zhang .  flag - flag indicating information about the preconditioner matrix
495aa644b49SHong Zhang           structure (same as flag in KSPSetOperators())
496aa644b49SHong Zhang -  ctx - [optional] user-defined context for matrix evaluation routine
497aa644b49SHong Zhang 
498aa644b49SHong Zhang    Notes:
499aa644b49SHong Zhang    See KSPSetOperators() for important information about setting the flag
500aa644b49SHong Zhang    output parameter in the routine func().  Be sure to read this information!
501aa644b49SHong Zhang 
502aa644b49SHong Zhang    The routine func() takes Mat * as the matrix arguments rather than Mat.
503aa644b49SHong Zhang    This allows the matrix evaluation routine to replace A and/or B with a
504aa644b49SHong Zhang    completely new new matrix structure (not just different matrix elements)
505aa644b49SHong Zhang    when appropriate, for instance, if the nonzero structure is changing
506aa644b49SHong Zhang    throughout the global iterations.
507aa644b49SHong Zhang 
508aa644b49SHong Zhang    Notes:
509aa644b49SHong Zhang    Currently, TSSetLHSMatrix() only supports the TS_BEULER type.
510aa644b49SHong Zhang 
511aa644b49SHong Zhang    Level: beginner
512aa644b49SHong Zhang 
513aa644b49SHong Zhang .keywords: TS, timestep, set, left-hand-side, matrix
514aa644b49SHong Zhang 
515aa644b49SHong Zhang .seealso: TSSetRHSMatrix()
516aa644b49SHong Zhang @*/
517aa644b49SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSetLHSMatrix(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
518aa644b49SHong Zhang {
519aa644b49SHong Zhang   PetscTruth     sametype;
520aa644b49SHong Zhang   PetscErrorCode ierr;
521aa644b49SHong Zhang 
522aa644b49SHong Zhang   PetscFunctionBegin;
523aa644b49SHong Zhang   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
524aa644b49SHong Zhang   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
525aa644b49SHong Zhang   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
526aa644b49SHong Zhang   PetscCheckSameComm(ts,1,A,2);
527aa644b49SHong Zhang   PetscCheckSameComm(ts,1,B,3);
528aa644b49SHong Zhang 
529aa644b49SHong Zhang   if (!ts->type_name) SETERRQ(PETSC_ERR_ARG_NULL,"TS type must be set before calling TSSetLHSMatrix()");
530aa644b49SHong Zhang   ierr = PetscTypeCompare((PetscObject)ts,"beuler",&sametype);CHKERRQ(ierr);
531aa644b49SHong Zhang   if (!sametype)
532aa644b49SHong Zhang     SETERRQ1(PETSC_ERR_SUP,"TS type %s not supported for LHSMatrix",ts->type_name);
533aa644b49SHong Zhang   if (ts->problem_type == TS_NONLINEAR) {
534aa644b49SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems yet");
535aa644b49SHong Zhang   }
536aa644b49SHong Zhang 
537aa644b49SHong Zhang   ts->ops->lhsmatrix = f;
538aa644b49SHong Zhang   ts->jacPlhs        = ctx;
539aa644b49SHong Zhang   ts->Alhs           = A;
540aa644b49SHong Zhang   PetscFunctionReturn(0);
541aa644b49SHong Zhang }
542aa644b49SHong Zhang 
543aa644b49SHong Zhang #undef __FUNCT__
5444a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
545d763cef2SBarry Smith /*@C
546d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
547d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
5483c94ec11SBarry Smith    Use TSSetRHSMatrix() for linear problems.
549d763cef2SBarry Smith 
550d763cef2SBarry Smith    Collective on TS
551d763cef2SBarry Smith 
552d763cef2SBarry Smith    Input Parameters:
553d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
554d763cef2SBarry Smith .  A   - Jacobian matrix
555d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
556d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
557d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
558d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
559d763cef2SBarry Smith 
560d763cef2SBarry Smith    Calling sequence of func:
56187828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
562d763cef2SBarry Smith 
563d763cef2SBarry Smith +  t - current timestep
564d763cef2SBarry Smith .  u - input vector
565d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
566d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
567d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
56894b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
569d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
570d763cef2SBarry Smith 
571d763cef2SBarry Smith    Notes:
57294b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
573d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
574d763cef2SBarry Smith 
575d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
576d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
577d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
578d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
579d763cef2SBarry Smith    throughout the global iterations.
580d763cef2SBarry Smith 
581d763cef2SBarry Smith    Level: beginner
582d763cef2SBarry Smith 
583d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
584d763cef2SBarry Smith 
585d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
5863c94ec11SBarry Smith           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetRHSMatrix()
587d763cef2SBarry Smith 
588d763cef2SBarry Smith @*/
58963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
590d763cef2SBarry Smith {
591d763cef2SBarry Smith   PetscFunctionBegin;
5924482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
5934482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
5944482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
595c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
596c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
597d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
59829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
599d763cef2SBarry Smith   }
600d763cef2SBarry Smith 
601000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
602d763cef2SBarry Smith   ts->jacP             = ctx;
6038beb423aSHong Zhang   ts->Arhs             = A;
604d763cef2SBarry Smith   ts->B                = B;
605d763cef2SBarry Smith   PetscFunctionReturn(0);
606d763cef2SBarry Smith }
607d763cef2SBarry Smith 
6084a2ae208SSatish Balay #undef __FUNCT__
6094a2ae208SSatish Balay #define __FUNCT__ "TSView"
6107e2c5f70SBarry Smith /*@C
611d763cef2SBarry Smith     TSView - Prints the TS data structure.
612d763cef2SBarry Smith 
6134c49b128SBarry Smith     Collective on TS
614d763cef2SBarry Smith 
615d763cef2SBarry Smith     Input Parameters:
616d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
617d763cef2SBarry Smith -   viewer - visualization context
618d763cef2SBarry Smith 
619d763cef2SBarry Smith     Options Database Key:
620d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
621d763cef2SBarry Smith 
622d763cef2SBarry Smith     Notes:
623d763cef2SBarry Smith     The available visualization contexts include
624b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
625b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
626d763cef2SBarry Smith          output where only the first processor opens
627d763cef2SBarry Smith          the file.  All other processors send their
628d763cef2SBarry Smith          data to the first processor to print.
629d763cef2SBarry Smith 
630d763cef2SBarry Smith     The user can open an alternative visualization context with
631b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
632d763cef2SBarry Smith 
633d763cef2SBarry Smith     Level: beginner
634d763cef2SBarry Smith 
635d763cef2SBarry Smith .keywords: TS, timestep, view
636d763cef2SBarry Smith 
637b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
638d763cef2SBarry Smith @*/
63963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer)
640d763cef2SBarry Smith {
641dfbe8321SBarry Smith   PetscErrorCode ierr;
642454a90a3SBarry Smith   char           *type;
64332077d6dSBarry Smith   PetscTruth     iascii,isstring;
644d763cef2SBarry Smith 
645d763cef2SBarry Smith   PetscFunctionBegin;
6464482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
6473050cee2SBarry Smith   if (!viewer) {
6483050cee2SBarry Smith     ierr = PetscViewerASCIIGetStdout(ts->comm,&viewer);CHKERRQ(ierr);
6493050cee2SBarry Smith   }
6504482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
651c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
652fd16b177SBarry Smith 
65332077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
654b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
65532077d6dSBarry Smith   if (iascii) {
656b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
657454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
658454a90a3SBarry Smith     if (type) {
659b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
660184914b5SBarry Smith     } else {
661b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
662184914b5SBarry Smith     }
663000e7ae3SMatthew Knepley     if (ts->ops->view) {
664b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
665000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
666b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
667d763cef2SBarry Smith     }
66877431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
669a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
670d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
67177431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
672d763cef2SBarry Smith     }
67377431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
6740f5bd95cSBarry Smith   } else if (isstring) {
675454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
676b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
677d763cef2SBarry Smith   }
678b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
67994b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
680d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
681b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
682d763cef2SBarry Smith   PetscFunctionReturn(0);
683d763cef2SBarry Smith }
684d763cef2SBarry Smith 
685d763cef2SBarry Smith 
6864a2ae208SSatish Balay #undef __FUNCT__
6874a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
688d763cef2SBarry Smith /*@C
689d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
690d763cef2SBarry Smith    the timesteppers.
691d763cef2SBarry Smith 
692d763cef2SBarry Smith    Collective on TS
693d763cef2SBarry Smith 
694d763cef2SBarry Smith    Input Parameters:
695d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
696d763cef2SBarry Smith -  usrP - optional user context
697d763cef2SBarry Smith 
698d763cef2SBarry Smith    Level: intermediate
699d763cef2SBarry Smith 
700d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
701d763cef2SBarry Smith 
702d763cef2SBarry Smith .seealso: TSGetApplicationContext()
703d763cef2SBarry Smith @*/
70463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP)
705d763cef2SBarry Smith {
706d763cef2SBarry Smith   PetscFunctionBegin;
7074482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
708d763cef2SBarry Smith   ts->user = usrP;
709d763cef2SBarry Smith   PetscFunctionReturn(0);
710d763cef2SBarry Smith }
711d763cef2SBarry Smith 
7124a2ae208SSatish Balay #undef __FUNCT__
7134a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
714d763cef2SBarry Smith /*@C
715d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
716d763cef2SBarry Smith     timestepper.
717d763cef2SBarry Smith 
718d763cef2SBarry Smith     Not Collective
719d763cef2SBarry Smith 
720d763cef2SBarry Smith     Input Parameter:
721d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
722d763cef2SBarry Smith 
723d763cef2SBarry Smith     Output Parameter:
724d763cef2SBarry Smith .   usrP - user context
725d763cef2SBarry Smith 
726d763cef2SBarry Smith     Level: intermediate
727d763cef2SBarry Smith 
728d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
729d763cef2SBarry Smith 
730d763cef2SBarry Smith .seealso: TSSetApplicationContext()
731d763cef2SBarry Smith @*/
73263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP)
733d763cef2SBarry Smith {
734d763cef2SBarry Smith   PetscFunctionBegin;
7354482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
736d763cef2SBarry Smith   *usrP = ts->user;
737d763cef2SBarry Smith   PetscFunctionReturn(0);
738d763cef2SBarry Smith }
739d763cef2SBarry Smith 
7404a2ae208SSatish Balay #undef __FUNCT__
7414a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
742d763cef2SBarry Smith /*@
743d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
744d763cef2SBarry Smith 
745d763cef2SBarry Smith    Not Collective
746d763cef2SBarry Smith 
747d763cef2SBarry Smith    Input Parameter:
748d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
749d763cef2SBarry Smith 
750d763cef2SBarry Smith    Output Parameter:
751d763cef2SBarry Smith .  iter - number steps so far
752d763cef2SBarry Smith 
753d763cef2SBarry Smith    Level: intermediate
754d763cef2SBarry Smith 
755d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
756d763cef2SBarry Smith @*/
75763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter)
758d763cef2SBarry Smith {
759d763cef2SBarry Smith   PetscFunctionBegin;
7604482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
7614482741eSBarry Smith   PetscValidIntPointer(iter,2);
762d763cef2SBarry Smith   *iter = ts->steps;
763d763cef2SBarry Smith   PetscFunctionReturn(0);
764d763cef2SBarry Smith }
765d763cef2SBarry Smith 
7664a2ae208SSatish Balay #undef __FUNCT__
7674a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
768d763cef2SBarry Smith /*@
769d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
770d763cef2SBarry Smith    as well as the initial time.
771d763cef2SBarry Smith 
772d763cef2SBarry Smith    Collective on TS
773d763cef2SBarry Smith 
774d763cef2SBarry Smith    Input Parameters:
775d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
776d763cef2SBarry Smith .  initial_time - the initial time
777d763cef2SBarry Smith -  time_step - the size of the timestep
778d763cef2SBarry Smith 
779d763cef2SBarry Smith    Level: intermediate
780d763cef2SBarry Smith 
781d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
782d763cef2SBarry Smith 
783d763cef2SBarry Smith .keywords: TS, set, initial, timestep
784d763cef2SBarry Smith @*/
78563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
786d763cef2SBarry Smith {
787d763cef2SBarry Smith   PetscFunctionBegin;
7884482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
789d763cef2SBarry Smith   ts->time_step         = time_step;
790d763cef2SBarry Smith   ts->initial_time_step = time_step;
791d763cef2SBarry Smith   ts->ptime             = initial_time;
792d763cef2SBarry Smith   PetscFunctionReturn(0);
793d763cef2SBarry Smith }
794d763cef2SBarry Smith 
7954a2ae208SSatish Balay #undef __FUNCT__
7964a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
797d763cef2SBarry Smith /*@
798d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
799d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
800d763cef2SBarry Smith 
801d763cef2SBarry Smith    Collective on TS
802d763cef2SBarry Smith 
803d763cef2SBarry Smith    Input Parameters:
804d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
805d763cef2SBarry Smith -  time_step - the size of the timestep
806d763cef2SBarry Smith 
807d763cef2SBarry Smith    Level: intermediate
808d763cef2SBarry Smith 
809d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
810d763cef2SBarry Smith 
811d763cef2SBarry Smith .keywords: TS, set, timestep
812d763cef2SBarry Smith @*/
81363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step)
814d763cef2SBarry Smith {
815d763cef2SBarry Smith   PetscFunctionBegin;
8164482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
817d763cef2SBarry Smith   ts->time_step = time_step;
818d763cef2SBarry Smith   PetscFunctionReturn(0);
819d763cef2SBarry Smith }
820d763cef2SBarry Smith 
8214a2ae208SSatish Balay #undef __FUNCT__
8224a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
823d763cef2SBarry Smith /*@
824d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
825d763cef2SBarry Smith 
826d763cef2SBarry Smith    Not Collective
827d763cef2SBarry Smith 
828d763cef2SBarry Smith    Input Parameter:
829d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
830d763cef2SBarry Smith 
831d763cef2SBarry Smith    Output Parameter:
832d763cef2SBarry Smith .  dt - the current timestep size
833d763cef2SBarry Smith 
834d763cef2SBarry Smith    Level: intermediate
835d763cef2SBarry Smith 
836d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
837d763cef2SBarry Smith 
838d763cef2SBarry Smith .keywords: TS, get, timestep
839d763cef2SBarry Smith @*/
84063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt)
841d763cef2SBarry Smith {
842d763cef2SBarry Smith   PetscFunctionBegin;
8434482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
8444482741eSBarry Smith   PetscValidDoublePointer(dt,2);
845d763cef2SBarry Smith   *dt = ts->time_step;
846d763cef2SBarry Smith   PetscFunctionReturn(0);
847d763cef2SBarry Smith }
848d763cef2SBarry Smith 
8494a2ae208SSatish Balay #undef __FUNCT__
8504a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
851d8e5e3e6SSatish Balay /*@
852d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
853d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
854d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
855d763cef2SBarry Smith    the solution at the next timestep has been calculated.
856d763cef2SBarry Smith 
857d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
858d763cef2SBarry Smith 
859d763cef2SBarry Smith    Input Parameter:
860d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
861d763cef2SBarry Smith 
862d763cef2SBarry Smith    Output Parameter:
863d763cef2SBarry Smith .  v - the vector containing the solution
864d763cef2SBarry Smith 
865d763cef2SBarry Smith    Level: intermediate
866d763cef2SBarry Smith 
867d763cef2SBarry Smith .seealso: TSGetTimeStep()
868d763cef2SBarry Smith 
869d763cef2SBarry Smith .keywords: TS, timestep, get, solution
870d763cef2SBarry Smith @*/
87163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v)
872d763cef2SBarry Smith {
873d763cef2SBarry Smith   PetscFunctionBegin;
8744482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
8754482741eSBarry Smith   PetscValidPointer(v,2);
876d763cef2SBarry Smith   *v = ts->vec_sol_always;
877d763cef2SBarry Smith   PetscFunctionReturn(0);
878d763cef2SBarry Smith }
879d763cef2SBarry Smith 
880bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
8814a2ae208SSatish Balay #undef __FUNCT__
882bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
883d8e5e3e6SSatish Balay /*@
884bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
885d763cef2SBarry Smith 
886bdad233fSMatthew Knepley   Not collective
887d763cef2SBarry Smith 
888bdad233fSMatthew Knepley   Input Parameters:
889bdad233fSMatthew Knepley + ts   - The TS
890bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
891d763cef2SBarry Smith .vb
892d763cef2SBarry Smith          U_t = A U
893d763cef2SBarry Smith          U_t = A(t) U
894d763cef2SBarry Smith          U_t = F(t,U)
895d763cef2SBarry Smith .ve
896d763cef2SBarry Smith 
897d763cef2SBarry Smith    Level: beginner
898d763cef2SBarry Smith 
899bdad233fSMatthew Knepley .keywords: TS, problem type
900bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
901d763cef2SBarry Smith @*/
90263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type)
903a7cc72afSBarry Smith {
904d763cef2SBarry Smith   PetscFunctionBegin;
9054482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
906bdad233fSMatthew Knepley   ts->problem_type = type;
907d763cef2SBarry Smith   PetscFunctionReturn(0);
908d763cef2SBarry Smith }
909d763cef2SBarry Smith 
910bdad233fSMatthew Knepley #undef __FUNCT__
911bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
912bdad233fSMatthew Knepley /*@C
913bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
914bdad233fSMatthew Knepley 
915bdad233fSMatthew Knepley   Not collective
916bdad233fSMatthew Knepley 
917bdad233fSMatthew Knepley   Input Parameter:
918bdad233fSMatthew Knepley . ts   - The TS
919bdad233fSMatthew Knepley 
920bdad233fSMatthew Knepley   Output Parameter:
921bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
922bdad233fSMatthew Knepley .vb
923bdad233fSMatthew Knepley          U_t = A U
924bdad233fSMatthew Knepley          U_t = A(t) U
925bdad233fSMatthew Knepley          U_t = F(t,U)
926bdad233fSMatthew Knepley .ve
927bdad233fSMatthew Knepley 
928bdad233fSMatthew Knepley    Level: beginner
929bdad233fSMatthew Knepley 
930bdad233fSMatthew Knepley .keywords: TS, problem type
931bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
932bdad233fSMatthew Knepley @*/
93363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type)
934a7cc72afSBarry Smith {
935bdad233fSMatthew Knepley   PetscFunctionBegin;
9364482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
9374482741eSBarry Smith   PetscValidIntPointer(type,2);
938bdad233fSMatthew Knepley   *type = ts->problem_type;
939bdad233fSMatthew Knepley   PetscFunctionReturn(0);
940bdad233fSMatthew Knepley }
941d763cef2SBarry Smith 
9424a2ae208SSatish Balay #undef __FUNCT__
9434a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
944d763cef2SBarry Smith /*@
945d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
946d763cef2SBarry Smith    of a timestepper.
947d763cef2SBarry Smith 
948d763cef2SBarry Smith    Collective on TS
949d763cef2SBarry Smith 
950d763cef2SBarry Smith    Input Parameter:
951d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
952d763cef2SBarry Smith 
953d763cef2SBarry Smith    Notes:
954d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
955d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
956d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
957d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
958d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
959d763cef2SBarry Smith 
960d763cef2SBarry Smith    Level: advanced
961d763cef2SBarry Smith 
962d763cef2SBarry Smith .keywords: TS, timestep, setup
963d763cef2SBarry Smith 
964d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
965d763cef2SBarry Smith @*/
96663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts)
967d763cef2SBarry Smith {
968dfbe8321SBarry Smith   PetscErrorCode ierr;
969d763cef2SBarry Smith 
970d763cef2SBarry Smith   PetscFunctionBegin;
9714482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
97229bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
973d763cef2SBarry Smith   if (!ts->type_name) {
974d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
975d763cef2SBarry Smith   }
976000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
977d763cef2SBarry Smith   ts->setupcalled = 1;
978d763cef2SBarry Smith   PetscFunctionReturn(0);
979d763cef2SBarry Smith }
980d763cef2SBarry Smith 
9814a2ae208SSatish Balay #undef __FUNCT__
9824a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
983d8e5e3e6SSatish Balay /*@
984d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
985d763cef2SBarry Smith    with TSCreate().
986d763cef2SBarry Smith 
987d763cef2SBarry Smith    Collective on TS
988d763cef2SBarry Smith 
989d763cef2SBarry Smith    Input Parameter:
990d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
991d763cef2SBarry Smith 
992d763cef2SBarry Smith    Level: beginner
993d763cef2SBarry Smith 
994d763cef2SBarry Smith .keywords: TS, timestepper, destroy
995d763cef2SBarry Smith 
996d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
997d763cef2SBarry Smith @*/
99863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts)
999d763cef2SBarry Smith {
10006849ba73SBarry Smith   PetscErrorCode ierr;
1001d763cef2SBarry Smith 
1002d763cef2SBarry Smith   PetscFunctionBegin;
10034482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1004d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
1005d763cef2SBarry Smith 
1006be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
10070f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
10088beb423aSHong Zhang   if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr)}
100994b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
1010d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
10111e3347e8SBarry Smith   if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);}
1012a6570f20SBarry Smith   ierr = TSMonitorCancel(ts);CHKERRQ(ierr);
1013a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1014d763cef2SBarry Smith   PetscFunctionReturn(0);
1015d763cef2SBarry Smith }
1016d763cef2SBarry Smith 
10174a2ae208SSatish Balay #undef __FUNCT__
10184a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
1019d8e5e3e6SSatish Balay /*@
1020d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1021d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
1022d763cef2SBarry Smith 
1023d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
1024d763cef2SBarry Smith 
1025d763cef2SBarry Smith    Input Parameter:
1026d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1027d763cef2SBarry Smith 
1028d763cef2SBarry Smith    Output Parameter:
1029d763cef2SBarry Smith .  snes - the nonlinear solver context
1030d763cef2SBarry Smith 
1031d763cef2SBarry Smith    Notes:
1032d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
1033d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
103494b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
1035d763cef2SBarry Smith 
1036d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
1037d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
1038d763cef2SBarry Smith 
1039d763cef2SBarry Smith    Level: beginner
1040d763cef2SBarry Smith 
1041d763cef2SBarry Smith .keywords: timestep, get, SNES
1042d763cef2SBarry Smith @*/
104363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes)
1044d763cef2SBarry Smith {
1045d763cef2SBarry Smith   PetscFunctionBegin;
10464482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
10474482741eSBarry Smith   PetscValidPointer(snes,2);
104894b7f48cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
1049d763cef2SBarry Smith   *snes = ts->snes;
1050d763cef2SBarry Smith   PetscFunctionReturn(0);
1051d763cef2SBarry Smith }
1052d763cef2SBarry Smith 
10534a2ae208SSatish Balay #undef __FUNCT__
105494b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
1055d8e5e3e6SSatish Balay /*@
105694b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
1057d763cef2SBarry Smith    a TS (timestepper) context.
1058d763cef2SBarry Smith 
105994b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1060d763cef2SBarry Smith 
1061d763cef2SBarry Smith    Input Parameter:
1062d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1063d763cef2SBarry Smith 
1064d763cef2SBarry Smith    Output Parameter:
106594b7f48cSBarry Smith .  ksp - the nonlinear solver context
1066d763cef2SBarry Smith 
1067d763cef2SBarry Smith    Notes:
106894b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1069d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1070d763cef2SBarry Smith    KSP and PC contexts as well.
1071d763cef2SBarry Smith 
107294b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
107394b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1074d763cef2SBarry Smith 
1075d763cef2SBarry Smith    Level: beginner
1076d763cef2SBarry Smith 
107794b7f48cSBarry Smith .keywords: timestep, get, KSP
1078d763cef2SBarry Smith @*/
107963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp)
1080d763cef2SBarry Smith {
1081d763cef2SBarry Smith   PetscFunctionBegin;
10824482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
10834482741eSBarry Smith   PetscValidPointer(ksp,2);
108429bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
108594b7f48cSBarry Smith   *ksp = ts->ksp;
1086d763cef2SBarry Smith   PetscFunctionReturn(0);
1087d763cef2SBarry Smith }
1088d763cef2SBarry Smith 
1089d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1090d763cef2SBarry Smith 
10914a2ae208SSatish Balay #undef __FUNCT__
1092adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1093adb62b0dSMatthew Knepley /*@
1094adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1095adb62b0dSMatthew Knepley    maximum time for iteration.
1096adb62b0dSMatthew Knepley 
1097adb62b0dSMatthew Knepley    Collective on TS
1098adb62b0dSMatthew Knepley 
1099adb62b0dSMatthew Knepley    Input Parameters:
1100adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1101adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1102adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1103adb62b0dSMatthew Knepley 
1104adb62b0dSMatthew Knepley    Level: intermediate
1105adb62b0dSMatthew Knepley 
1106adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1107adb62b0dSMatthew Knepley @*/
110863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1109adb62b0dSMatthew Knepley {
1110adb62b0dSMatthew Knepley   PetscFunctionBegin;
11114482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1112abc0a331SBarry Smith   if (maxsteps) {
11134482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1114adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1115adb62b0dSMatthew Knepley   }
1116abc0a331SBarry Smith   if (maxtime ) {
11174482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1118adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1119adb62b0dSMatthew Knepley   }
1120adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1121adb62b0dSMatthew Knepley }
1122adb62b0dSMatthew Knepley 
1123adb62b0dSMatthew Knepley #undef __FUNCT__
11244a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1125d763cef2SBarry Smith /*@
1126d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1127d763cef2SBarry Smith    maximum time for iteration.
1128d763cef2SBarry Smith 
1129d763cef2SBarry Smith    Collective on TS
1130d763cef2SBarry Smith 
1131d763cef2SBarry Smith    Input Parameters:
1132d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1133d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1134d763cef2SBarry Smith -  maxtime - final time to iterate to
1135d763cef2SBarry Smith 
1136d763cef2SBarry Smith    Options Database Keys:
1137d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1138d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1139d763cef2SBarry Smith 
1140d763cef2SBarry Smith    Notes:
1141d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1142d763cef2SBarry Smith 
1143d763cef2SBarry Smith    Level: intermediate
1144d763cef2SBarry Smith 
1145d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1146d763cef2SBarry Smith @*/
114763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1148d763cef2SBarry Smith {
1149d763cef2SBarry Smith   PetscFunctionBegin;
11504482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1151d763cef2SBarry Smith   ts->max_steps = maxsteps;
1152d763cef2SBarry Smith   ts->max_time  = maxtime;
1153d763cef2SBarry Smith   PetscFunctionReturn(0);
1154d763cef2SBarry Smith }
1155d763cef2SBarry Smith 
11564a2ae208SSatish Balay #undef __FUNCT__
11574a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1158d763cef2SBarry Smith /*@
1159d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1160d763cef2SBarry Smith    for use by the TS routines.
1161d763cef2SBarry Smith 
1162d763cef2SBarry Smith    Collective on TS and Vec
1163d763cef2SBarry Smith 
1164d763cef2SBarry Smith    Input Parameters:
1165d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1166d763cef2SBarry Smith -  x - the solution vector
1167d763cef2SBarry Smith 
1168d763cef2SBarry Smith    Level: beginner
1169d763cef2SBarry Smith 
1170d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1171d763cef2SBarry Smith @*/
117263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x)
1173d763cef2SBarry Smith {
1174d763cef2SBarry Smith   PetscFunctionBegin;
11754482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
11764482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1177d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1178d763cef2SBarry Smith   PetscFunctionReturn(0);
1179d763cef2SBarry Smith }
1180d763cef2SBarry Smith 
1181e74ef692SMatthew Knepley #undef __FUNCT__
1182e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1183ac226902SBarry Smith /*@C
1184000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1185000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1186000e7ae3SMatthew Knepley 
1187000e7ae3SMatthew Knepley   Collective on TS
1188000e7ae3SMatthew Knepley 
1189000e7ae3SMatthew Knepley   Input Parameters:
1190000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1191000e7ae3SMatthew Knepley - func - The function
1192000e7ae3SMatthew Knepley 
1193000e7ae3SMatthew Knepley   Calling sequence of func:
1194000e7ae3SMatthew Knepley . func (TS ts);
1195000e7ae3SMatthew Knepley 
1196000e7ae3SMatthew Knepley   Level: intermediate
1197000e7ae3SMatthew Knepley 
1198000e7ae3SMatthew Knepley .keywords: TS, timestep
1199000e7ae3SMatthew Knepley @*/
120063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1201000e7ae3SMatthew Knepley {
1202000e7ae3SMatthew Knepley   PetscFunctionBegin;
12034482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1204000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1205000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1206000e7ae3SMatthew Knepley }
1207000e7ae3SMatthew Knepley 
1208e74ef692SMatthew Knepley #undef __FUNCT__
1209e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1210000e7ae3SMatthew Knepley /*@
1211000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1212000e7ae3SMatthew Knepley 
1213000e7ae3SMatthew Knepley   Collective on TS
1214000e7ae3SMatthew Knepley 
1215000e7ae3SMatthew Knepley   Input Parameters:
1216000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1217000e7ae3SMatthew Knepley 
1218000e7ae3SMatthew Knepley   Level: developer
1219000e7ae3SMatthew Knepley 
1220000e7ae3SMatthew Knepley .keywords: TS, timestep
1221000e7ae3SMatthew Knepley @*/
122263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts)
1223000e7ae3SMatthew Knepley {
1224000e7ae3SMatthew Knepley   PetscFunctionBegin;
1225000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1226000e7ae3SMatthew Knepley }
1227000e7ae3SMatthew Knepley 
1228e74ef692SMatthew Knepley #undef __FUNCT__
1229e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1230ac226902SBarry Smith /*@C
1231000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1232000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1233000e7ae3SMatthew Knepley   the time step.
1234000e7ae3SMatthew Knepley 
1235000e7ae3SMatthew Knepley   Collective on TS
1236000e7ae3SMatthew Knepley 
1237000e7ae3SMatthew Knepley   Input Parameters:
1238000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1239000e7ae3SMatthew Knepley - func - The function
1240000e7ae3SMatthew Knepley 
1241000e7ae3SMatthew Knepley   Calling sequence of func:
1242000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1243000e7ae3SMatthew Knepley 
1244000e7ae3SMatthew Knepley + t   - The current time
1245000e7ae3SMatthew Knepley - dt  - The current time step
1246000e7ae3SMatthew Knepley 
1247000e7ae3SMatthew Knepley   Level: intermediate
1248000e7ae3SMatthew Knepley 
1249000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1250000e7ae3SMatthew Knepley @*/
125163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUpdate(TS ts, PetscErrorCode (*func)(TS, PetscReal, PetscReal *))
1252000e7ae3SMatthew Knepley {
1253000e7ae3SMatthew Knepley   PetscFunctionBegin;
12544482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1255000e7ae3SMatthew Knepley   ts->ops->update = func;
1256000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1257000e7ae3SMatthew Knepley }
1258000e7ae3SMatthew Knepley 
1259e74ef692SMatthew Knepley #undef __FUNCT__
1260e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1261000e7ae3SMatthew Knepley /*@
1262000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1263000e7ae3SMatthew Knepley 
1264000e7ae3SMatthew Knepley   Collective on TS
1265000e7ae3SMatthew Knepley 
1266000e7ae3SMatthew Knepley   Input Parameters:
1267000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1268000e7ae3SMatthew Knepley - t   - The current time
1269000e7ae3SMatthew Knepley 
1270000e7ae3SMatthew Knepley   Output Parameters:
1271000e7ae3SMatthew Knepley . dt  - The current time step
1272000e7ae3SMatthew Knepley 
1273000e7ae3SMatthew Knepley   Level: developer
1274000e7ae3SMatthew Knepley 
1275000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1276000e7ae3SMatthew Knepley @*/
127763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultUpdate(TS ts, PetscReal t, PetscReal *dt)
1278000e7ae3SMatthew Knepley {
1279000e7ae3SMatthew Knepley   PetscFunctionBegin;
1280000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1281000e7ae3SMatthew Knepley }
1282000e7ae3SMatthew Knepley 
1283e74ef692SMatthew Knepley #undef __FUNCT__
1284e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1285ac226902SBarry Smith /*@C
1286000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1287000e7ae3SMatthew Knepley   called once at the end of time stepping.
1288000e7ae3SMatthew Knepley 
1289000e7ae3SMatthew Knepley   Collective on TS
1290000e7ae3SMatthew Knepley 
1291000e7ae3SMatthew Knepley   Input Parameters:
1292000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1293000e7ae3SMatthew Knepley - func - The function
1294000e7ae3SMatthew Knepley 
1295000e7ae3SMatthew Knepley   Calling sequence of func:
1296000e7ae3SMatthew Knepley . func (TS ts);
1297000e7ae3SMatthew Knepley 
1298000e7ae3SMatthew Knepley   Level: intermediate
1299000e7ae3SMatthew Knepley 
1300000e7ae3SMatthew Knepley .keywords: TS, timestep
1301000e7ae3SMatthew Knepley @*/
130263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1303000e7ae3SMatthew Knepley {
1304000e7ae3SMatthew Knepley   PetscFunctionBegin;
13054482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1306000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1307000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1308000e7ae3SMatthew Knepley }
1309000e7ae3SMatthew Knepley 
1310e74ef692SMatthew Knepley #undef __FUNCT__
1311e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1312000e7ae3SMatthew Knepley /*@
1313000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1314000e7ae3SMatthew Knepley 
1315000e7ae3SMatthew Knepley   Collective on TS
1316000e7ae3SMatthew Knepley 
1317000e7ae3SMatthew Knepley   Input Parameters:
1318000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1319000e7ae3SMatthew Knepley 
1320000e7ae3SMatthew Knepley   Level: developer
1321000e7ae3SMatthew Knepley 
1322000e7ae3SMatthew Knepley .keywords: TS, timestep
1323000e7ae3SMatthew Knepley @*/
132463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts)
1325000e7ae3SMatthew Knepley {
1326000e7ae3SMatthew Knepley   PetscFunctionBegin;
1327000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1328000e7ae3SMatthew Knepley }
1329000e7ae3SMatthew Knepley 
1330d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1331d763cef2SBarry Smith 
13324a2ae208SSatish Balay #undef __FUNCT__
1333a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet"
1334d763cef2SBarry Smith /*@C
1335a6570f20SBarry Smith    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1336d763cef2SBarry Smith    timestep to display the iteration's  progress.
1337d763cef2SBarry Smith 
1338d763cef2SBarry Smith    Collective on TS
1339d763cef2SBarry Smith 
1340d763cef2SBarry Smith    Input Parameters:
1341d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1342d763cef2SBarry Smith .  func - monitoring routine
1343329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1344b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1345b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1346b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1347d763cef2SBarry Smith 
1348d763cef2SBarry Smith    Calling sequence of func:
1349a7cc72afSBarry Smith $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1350d763cef2SBarry Smith 
1351d763cef2SBarry Smith +    ts - the TS context
1352d763cef2SBarry Smith .    steps - iteration number
13531f06c33eSBarry Smith .    time - current time
1354d763cef2SBarry Smith .    x - current iterate
1355d763cef2SBarry Smith -    mctx - [optional] monitoring context
1356d763cef2SBarry Smith 
1357d763cef2SBarry Smith    Notes:
1358d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1359d763cef2SBarry Smith    already has been loaded.
1360d763cef2SBarry Smith 
1361d763cef2SBarry Smith    Level: intermediate
1362d763cef2SBarry Smith 
1363d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1364d763cef2SBarry Smith 
1365a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel()
1366d763cef2SBarry Smith @*/
1367a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1368d763cef2SBarry Smith {
1369d763cef2SBarry Smith   PetscFunctionBegin;
13704482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1371d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
137229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1373d763cef2SBarry Smith   }
1374d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1375329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1376d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1377d763cef2SBarry Smith   PetscFunctionReturn(0);
1378d763cef2SBarry Smith }
1379d763cef2SBarry Smith 
13804a2ae208SSatish Balay #undef __FUNCT__
1381a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel"
1382d763cef2SBarry Smith /*@C
1383a6570f20SBarry Smith    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1384d763cef2SBarry Smith 
1385d763cef2SBarry Smith    Collective on TS
1386d763cef2SBarry Smith 
1387d763cef2SBarry Smith    Input Parameters:
1388d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1389d763cef2SBarry Smith 
1390d763cef2SBarry Smith    Notes:
1391d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1392d763cef2SBarry Smith 
1393d763cef2SBarry Smith    Level: intermediate
1394d763cef2SBarry Smith 
1395d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1396d763cef2SBarry Smith 
1397a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet()
1398d763cef2SBarry Smith @*/
1399a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts)
1400d763cef2SBarry Smith {
1401d952e501SBarry Smith   PetscErrorCode ierr;
1402d952e501SBarry Smith   PetscInt       i;
1403d952e501SBarry Smith 
1404d763cef2SBarry Smith   PetscFunctionBegin;
14054482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1406d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1407d952e501SBarry Smith     if (ts->mdestroy[i]) {
1408d952e501SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
1409d952e501SBarry Smith     }
1410d952e501SBarry Smith   }
1411d763cef2SBarry Smith   ts->numbermonitors = 0;
1412d763cef2SBarry Smith   PetscFunctionReturn(0);
1413d763cef2SBarry Smith }
1414d763cef2SBarry Smith 
14154a2ae208SSatish Balay #undef __FUNCT__
1416a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault"
1417d8e5e3e6SSatish Balay /*@
1418a6570f20SBarry Smith    TSMonitorDefault - Sets the Default monitor
1419d8e5e3e6SSatish Balay @*/
1420a6570f20SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1421d763cef2SBarry Smith {
1422dfbe8321SBarry Smith   PetscErrorCode          ierr;
1423a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
1424d132466eSBarry Smith 
1425d763cef2SBarry Smith   PetscFunctionBegin;
1426a34d58ebSBarry Smith   if (!ctx) {
1427a34d58ebSBarry Smith     ierr = PetscViewerASCIIMonitorCreate(ts->comm,"stdout",0,&viewer);CHKERRQ(ierr);
1428a34d58ebSBarry Smith   }
1429a34d58ebSBarry Smith   ierr = PetscViewerASCIIMonitorPrintf(viewer,"timestep %D dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1430a34d58ebSBarry Smith   if (!ctx) {
1431a34d58ebSBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
1432a34d58ebSBarry Smith   }
1433d763cef2SBarry Smith   PetscFunctionReturn(0);
1434d763cef2SBarry Smith }
1435d763cef2SBarry Smith 
14364a2ae208SSatish Balay #undef __FUNCT__
14374a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1438d763cef2SBarry Smith /*@
1439d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1440d763cef2SBarry Smith 
1441d763cef2SBarry Smith    Collective on TS
1442d763cef2SBarry Smith 
1443d763cef2SBarry Smith    Input Parameter:
1444d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1445d763cef2SBarry Smith 
1446d763cef2SBarry Smith    Output Parameters:
1447d763cef2SBarry Smith +  steps - number of iterations until termination
1448142b95e3SSatish Balay -  ptime - time until termination
1449d763cef2SBarry Smith 
1450d763cef2SBarry Smith    Level: beginner
1451d763cef2SBarry Smith 
1452d763cef2SBarry Smith .keywords: TS, timestep, solve
1453d763cef2SBarry Smith 
1454d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1455d763cef2SBarry Smith @*/
145663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1457d763cef2SBarry Smith {
1458dfbe8321SBarry Smith   PetscErrorCode ierr;
1459d763cef2SBarry Smith 
1460d763cef2SBarry Smith   PetscFunctionBegin;
14614482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1462d405a339SMatthew Knepley   if (!ts->setupcalled) {
1463d405a339SMatthew Knepley     ierr = TSSetUp(ts);CHKERRQ(ierr);
1464d405a339SMatthew Knepley   }
1465d405a339SMatthew Knepley 
1466d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1467d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1468000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1469d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1470d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1471d405a339SMatthew Knepley 
14724bb05414SBarry Smith   if (!PetscPreLoadingOn) {
14734bb05414SBarry Smith     ierr = TSViewFromOptions(ts,ts->name);CHKERRQ(ierr);
1474d405a339SMatthew Knepley   }
1475d763cef2SBarry Smith   PetscFunctionReturn(0);
1476d763cef2SBarry Smith }
1477d763cef2SBarry Smith 
14784a2ae208SSatish Balay #undef __FUNCT__
14796a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve"
14806a4d4014SLisandro Dalcin /*@
14816a4d4014SLisandro Dalcin    TSSolve - Steps the requested number of timesteps.
14826a4d4014SLisandro Dalcin 
14836a4d4014SLisandro Dalcin    Collective on TS
14846a4d4014SLisandro Dalcin 
14856a4d4014SLisandro Dalcin    Input Parameter:
14866a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
14876a4d4014SLisandro Dalcin -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
14886a4d4014SLisandro Dalcin 
14896a4d4014SLisandro Dalcin    Level: beginner
14906a4d4014SLisandro Dalcin 
14916a4d4014SLisandro Dalcin .keywords: TS, timestep, solve
14926a4d4014SLisandro Dalcin 
14936a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep()
14946a4d4014SLisandro Dalcin @*/
14956a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x)
14966a4d4014SLisandro Dalcin {
14976a4d4014SLisandro Dalcin   PetscInt       steps;
14986a4d4014SLisandro Dalcin   PetscReal      ptime;
14996a4d4014SLisandro Dalcin   PetscErrorCode ierr;
15006a4d4014SLisandro Dalcin   PetscFunctionBegin;
15016a4d4014SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
15026a4d4014SLisandro Dalcin   /* set solution vector if provided */
15036a4d4014SLisandro Dalcin   if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); }
15046a4d4014SLisandro Dalcin   /* reset time step and iteration counters */
15056a4d4014SLisandro Dalcin   ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0;
15066a4d4014SLisandro Dalcin   /* steps the requested number of timesteps. */
15076a4d4014SLisandro Dalcin   ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr);
15086a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
15096a4d4014SLisandro Dalcin }
15106a4d4014SLisandro Dalcin 
15116a4d4014SLisandro Dalcin #undef __FUNCT__
15124a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1513d763cef2SBarry Smith /*
1514d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1515d763cef2SBarry Smith */
1516a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1517d763cef2SBarry Smith {
15186849ba73SBarry Smith   PetscErrorCode ierr;
1519a7cc72afSBarry Smith   PetscInt i,n = ts->numbermonitors;
1520d763cef2SBarry Smith 
1521d763cef2SBarry Smith   PetscFunctionBegin;
1522d763cef2SBarry Smith   for (i=0; i<n; i++) {
1523142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1524d763cef2SBarry Smith   }
1525d763cef2SBarry Smith   PetscFunctionReturn(0);
1526d763cef2SBarry Smith }
1527d763cef2SBarry Smith 
1528d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1529d763cef2SBarry Smith 
15304a2ae208SSatish Balay #undef __FUNCT__
1531a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate"
1532d763cef2SBarry Smith /*@C
1533a6570f20SBarry Smith    TSMonitorLGCreate - Creates a line graph context for use with
1534d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1535d763cef2SBarry Smith 
1536d763cef2SBarry Smith    Collective on TS
1537d763cef2SBarry Smith 
1538d763cef2SBarry Smith    Input Parameters:
1539d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1540d763cef2SBarry Smith .  label - the title to put in the title bar
15417c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1542d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1543d763cef2SBarry Smith 
1544d763cef2SBarry Smith    Output Parameter:
1545d763cef2SBarry Smith .  draw - the drawing context
1546d763cef2SBarry Smith 
1547d763cef2SBarry Smith    Options Database Key:
1548a6570f20SBarry Smith .  -ts_monitor_draw - automatically sets line graph monitor
1549d763cef2SBarry Smith 
1550d763cef2SBarry Smith    Notes:
1551a6570f20SBarry Smith    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1552d763cef2SBarry Smith 
1553d763cef2SBarry Smith    Level: intermediate
1554d763cef2SBarry Smith 
15557c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1556d763cef2SBarry Smith 
1557a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet()
15587c922b88SBarry Smith 
1559d763cef2SBarry Smith @*/
1560a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1561d763cef2SBarry Smith {
1562b0a32e0cSBarry Smith   PetscDraw      win;
1563dfbe8321SBarry Smith   PetscErrorCode ierr;
1564d763cef2SBarry Smith 
1565d763cef2SBarry Smith   PetscFunctionBegin;
1566b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1567b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1568b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1569b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1570d763cef2SBarry Smith 
157152e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1572d763cef2SBarry Smith   PetscFunctionReturn(0);
1573d763cef2SBarry Smith }
1574d763cef2SBarry Smith 
15754a2ae208SSatish Balay #undef __FUNCT__
1576a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG"
1577a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1578d763cef2SBarry Smith {
1579b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
158087828ca2SBarry Smith   PetscReal      x,y = ptime;
1581dfbe8321SBarry Smith   PetscErrorCode ierr;
1582d763cef2SBarry Smith 
1583d763cef2SBarry Smith   PetscFunctionBegin;
15847c922b88SBarry Smith   if (!monctx) {
15857c922b88SBarry Smith     MPI_Comm    comm;
1586b0a32e0cSBarry Smith     PetscViewer viewer;
15877c922b88SBarry Smith 
15887c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1589b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1590b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
15917c922b88SBarry Smith   }
15927c922b88SBarry Smith 
1593b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
159487828ca2SBarry Smith   x = (PetscReal)n;
1595b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1596d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1597b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1598d763cef2SBarry Smith   }
1599d763cef2SBarry Smith   PetscFunctionReturn(0);
1600d763cef2SBarry Smith }
1601d763cef2SBarry Smith 
16024a2ae208SSatish Balay #undef __FUNCT__
1603a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy"
1604d763cef2SBarry Smith /*@C
1605a6570f20SBarry Smith    TSMonitorLGDestroy - Destroys a line graph context that was created
1606a6570f20SBarry Smith    with TSMonitorLGCreate().
1607d763cef2SBarry Smith 
1608b0a32e0cSBarry Smith    Collective on PetscDrawLG
1609d763cef2SBarry Smith 
1610d763cef2SBarry Smith    Input Parameter:
1611d763cef2SBarry Smith .  draw - the drawing context
1612d763cef2SBarry Smith 
1613d763cef2SBarry Smith    Level: intermediate
1614d763cef2SBarry Smith 
1615d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1616d763cef2SBarry Smith 
1617a6570f20SBarry Smith .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1618d763cef2SBarry Smith @*/
1619a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg)
1620d763cef2SBarry Smith {
1621b0a32e0cSBarry Smith   PetscDraw      draw;
1622dfbe8321SBarry Smith   PetscErrorCode ierr;
1623d763cef2SBarry Smith 
1624d763cef2SBarry Smith   PetscFunctionBegin;
1625b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1626b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1627b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1628d763cef2SBarry Smith   PetscFunctionReturn(0);
1629d763cef2SBarry Smith }
1630d763cef2SBarry Smith 
16314a2ae208SSatish Balay #undef __FUNCT__
16324a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1633d763cef2SBarry Smith /*@
1634d763cef2SBarry Smith    TSGetTime - Gets the current time.
1635d763cef2SBarry Smith 
1636d763cef2SBarry Smith    Not Collective
1637d763cef2SBarry Smith 
1638d763cef2SBarry Smith    Input Parameter:
1639d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1640d763cef2SBarry Smith 
1641d763cef2SBarry Smith    Output Parameter:
1642d763cef2SBarry Smith .  t  - the current time
1643d763cef2SBarry Smith 
1644d763cef2SBarry Smith    Contributed by: Matthew Knepley
1645d763cef2SBarry Smith 
1646d763cef2SBarry Smith    Level: beginner
1647d763cef2SBarry Smith 
1648d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1649d763cef2SBarry Smith 
1650d763cef2SBarry Smith .keywords: TS, get, time
1651d763cef2SBarry Smith @*/
165263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t)
1653d763cef2SBarry Smith {
1654d763cef2SBarry Smith   PetscFunctionBegin;
16554482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
16564482741eSBarry Smith   PetscValidDoublePointer(t,2);
1657d763cef2SBarry Smith   *t = ts->ptime;
1658d763cef2SBarry Smith   PetscFunctionReturn(0);
1659d763cef2SBarry Smith }
1660d763cef2SBarry Smith 
16614a2ae208SSatish Balay #undef __FUNCT__
16626a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime"
16636a4d4014SLisandro Dalcin /*@
16646a4d4014SLisandro Dalcin    TSSetTime - Allows one to reset the time.
16656a4d4014SLisandro Dalcin 
16666a4d4014SLisandro Dalcin    Collective on TS
16676a4d4014SLisandro Dalcin 
16686a4d4014SLisandro Dalcin    Input Parameters:
16696a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
16706a4d4014SLisandro Dalcin -  time - the time
16716a4d4014SLisandro Dalcin 
16726a4d4014SLisandro Dalcin    Level: intermediate
16736a4d4014SLisandro Dalcin 
16746a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration()
16756a4d4014SLisandro Dalcin 
16766a4d4014SLisandro Dalcin .keywords: TS, set, time
16776a4d4014SLisandro Dalcin @*/
16786a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t)
16796a4d4014SLisandro Dalcin {
16806a4d4014SLisandro Dalcin   PetscFunctionBegin;
16816a4d4014SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
16826a4d4014SLisandro Dalcin   ts->ptime = t;
16836a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
16846a4d4014SLisandro Dalcin }
16856a4d4014SLisandro Dalcin 
16866a4d4014SLisandro Dalcin #undef __FUNCT__
16874a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1688d763cef2SBarry Smith /*@C
1689d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1690d763cef2SBarry Smith    TS options in the database.
1691d763cef2SBarry Smith 
1692d763cef2SBarry Smith    Collective on TS
1693d763cef2SBarry Smith 
1694d763cef2SBarry Smith    Input Parameter:
1695d763cef2SBarry Smith +  ts     - The TS context
1696d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1697d763cef2SBarry Smith 
1698d763cef2SBarry Smith    Notes:
1699d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1700d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1701d763cef2SBarry Smith    hyphen.
1702d763cef2SBarry Smith 
1703d763cef2SBarry Smith    Contributed by: Matthew Knepley
1704d763cef2SBarry Smith 
1705d763cef2SBarry Smith    Level: advanced
1706d763cef2SBarry Smith 
1707d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1708d763cef2SBarry Smith 
1709d763cef2SBarry Smith .seealso: TSSetFromOptions()
1710d763cef2SBarry Smith 
1711d763cef2SBarry Smith @*/
171263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[])
1713d763cef2SBarry Smith {
1714dfbe8321SBarry Smith   PetscErrorCode ierr;
1715d763cef2SBarry Smith 
1716d763cef2SBarry Smith   PetscFunctionBegin;
17174482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1718d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1719d763cef2SBarry Smith   switch(ts->problem_type) {
1720d763cef2SBarry Smith     case TS_NONLINEAR:
1721d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1722d763cef2SBarry Smith       break;
1723d763cef2SBarry Smith     case TS_LINEAR:
172494b7f48cSBarry Smith       ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1725d763cef2SBarry Smith       break;
1726d763cef2SBarry Smith   }
1727d763cef2SBarry Smith   PetscFunctionReturn(0);
1728d763cef2SBarry Smith }
1729d763cef2SBarry Smith 
1730d763cef2SBarry Smith 
17314a2ae208SSatish Balay #undef __FUNCT__
17324a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1733d763cef2SBarry Smith /*@C
1734d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1735d763cef2SBarry Smith    TS options in the database.
1736d763cef2SBarry Smith 
1737d763cef2SBarry Smith    Collective on TS
1738d763cef2SBarry Smith 
1739d763cef2SBarry Smith    Input Parameter:
1740d763cef2SBarry Smith +  ts     - The TS context
1741d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1742d763cef2SBarry Smith 
1743d763cef2SBarry Smith    Notes:
1744d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1745d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1746d763cef2SBarry Smith    hyphen.
1747d763cef2SBarry Smith 
1748d763cef2SBarry Smith    Contributed by: Matthew Knepley
1749d763cef2SBarry Smith 
1750d763cef2SBarry Smith    Level: advanced
1751d763cef2SBarry Smith 
1752d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1753d763cef2SBarry Smith 
1754d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1755d763cef2SBarry Smith 
1756d763cef2SBarry Smith @*/
175763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[])
1758d763cef2SBarry Smith {
1759dfbe8321SBarry Smith   PetscErrorCode ierr;
1760d763cef2SBarry Smith 
1761d763cef2SBarry Smith   PetscFunctionBegin;
17624482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1763d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1764d763cef2SBarry Smith   switch(ts->problem_type) {
1765d763cef2SBarry Smith     case TS_NONLINEAR:
1766d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1767d763cef2SBarry Smith       break;
1768d763cef2SBarry Smith     case TS_LINEAR:
176994b7f48cSBarry Smith       ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
1770d763cef2SBarry Smith       break;
1771d763cef2SBarry Smith   }
1772d763cef2SBarry Smith   PetscFunctionReturn(0);
1773d763cef2SBarry Smith }
1774d763cef2SBarry Smith 
17754a2ae208SSatish Balay #undef __FUNCT__
17764a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1777d763cef2SBarry Smith /*@C
1778d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1779d763cef2SBarry Smith    TS options in the database.
1780d763cef2SBarry Smith 
1781d763cef2SBarry Smith    Not Collective
1782d763cef2SBarry Smith 
1783d763cef2SBarry Smith    Input Parameter:
1784d763cef2SBarry Smith .  ts - The TS context
1785d763cef2SBarry Smith 
1786d763cef2SBarry Smith    Output Parameter:
1787d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1788d763cef2SBarry Smith 
1789d763cef2SBarry Smith    Contributed by: Matthew Knepley
1790d763cef2SBarry Smith 
1791d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1792d763cef2SBarry Smith    sufficient length to hold the prefix.
1793d763cef2SBarry Smith 
1794d763cef2SBarry Smith    Level: intermediate
1795d763cef2SBarry Smith 
1796d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1797d763cef2SBarry Smith 
1798d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1799d763cef2SBarry Smith @*/
1800e060cb09SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[])
1801d763cef2SBarry Smith {
1802dfbe8321SBarry Smith   PetscErrorCode ierr;
1803d763cef2SBarry Smith 
1804d763cef2SBarry Smith   PetscFunctionBegin;
18054482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
18064482741eSBarry Smith   PetscValidPointer(prefix,2);
1807d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1808d763cef2SBarry Smith   PetscFunctionReturn(0);
1809d763cef2SBarry Smith }
1810d763cef2SBarry Smith 
18114a2ae208SSatish Balay #undef __FUNCT__
18124a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1813d763cef2SBarry Smith /*@C
1814d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1815d763cef2SBarry Smith 
1816d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1817d763cef2SBarry Smith 
1818d763cef2SBarry Smith    Input Parameter:
1819d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1820d763cef2SBarry Smith 
1821d763cef2SBarry Smith    Output Parameters:
1822d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1823d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1824d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1825d763cef2SBarry Smith 
1826d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1827d763cef2SBarry Smith 
1828d763cef2SBarry Smith    Contributed by: Matthew Knepley
1829d763cef2SBarry Smith 
1830d763cef2SBarry Smith    Level: intermediate
1831d763cef2SBarry Smith 
1832d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1833d763cef2SBarry Smith 
1834d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1835d763cef2SBarry Smith 
1836d763cef2SBarry Smith @*/
183763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1838d763cef2SBarry Smith {
1839d763cef2SBarry Smith   PetscFunctionBegin;
18404482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
18418beb423aSHong Zhang   if (A)   *A = ts->Arhs;
1842d763cef2SBarry Smith   if (M)   *M = ts->B;
1843d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1844d763cef2SBarry Smith   PetscFunctionReturn(0);
1845d763cef2SBarry Smith }
1846d763cef2SBarry Smith 
18474a2ae208SSatish Balay #undef __FUNCT__
18484a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1849d763cef2SBarry Smith /*@C
1850d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1851d763cef2SBarry Smith 
1852d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1853d763cef2SBarry Smith 
1854d763cef2SBarry Smith    Input Parameter:
1855d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1856d763cef2SBarry Smith 
1857d763cef2SBarry Smith    Output Parameters:
1858d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1859d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1860d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1861d763cef2SBarry Smith 
1862d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1863d763cef2SBarry Smith 
1864d763cef2SBarry Smith    Contributed by: Matthew Knepley
1865d763cef2SBarry Smith 
1866d763cef2SBarry Smith    Level: intermediate
1867d763cef2SBarry Smith 
1868d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1869d763cef2SBarry Smith 
1870d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1871d763cef2SBarry Smith @*/
187263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1873d763cef2SBarry Smith {
1874dfbe8321SBarry Smith   PetscErrorCode ierr;
1875d763cef2SBarry Smith 
1876d763cef2SBarry Smith   PetscFunctionBegin;
1877d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1878d763cef2SBarry Smith   PetscFunctionReturn(0);
1879d763cef2SBarry Smith }
1880d763cef2SBarry Smith 
18811713a123SBarry Smith #undef __FUNCT__
1882a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution"
18831713a123SBarry Smith /*@C
1884a6570f20SBarry Smith    TSMonitorSolution - Monitors progress of the TS solvers by calling
18851713a123SBarry Smith    VecView() for the solution at each timestep
18861713a123SBarry Smith 
18871713a123SBarry Smith    Collective on TS
18881713a123SBarry Smith 
18891713a123SBarry Smith    Input Parameters:
18901713a123SBarry Smith +  ts - the TS context
18911713a123SBarry Smith .  step - current time-step
1892142b95e3SSatish Balay .  ptime - current time
18931713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
18941713a123SBarry Smith 
18951713a123SBarry Smith    Level: intermediate
18961713a123SBarry Smith 
18971713a123SBarry Smith .keywords: TS,  vector, monitor, view
18981713a123SBarry Smith 
1899a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
19001713a123SBarry Smith @*/
1901a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
19021713a123SBarry Smith {
1903dfbe8321SBarry Smith   PetscErrorCode ierr;
19041713a123SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
19051713a123SBarry Smith 
19061713a123SBarry Smith   PetscFunctionBegin;
1907a34d58ebSBarry Smith   if (!dummy) {
1908a34d58ebSBarry Smith     viewer = PETSC_VIEWER_DRAW_(ts->comm);
19091713a123SBarry Smith   }
19101713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
19111713a123SBarry Smith   PetscFunctionReturn(0);
19121713a123SBarry Smith }
19131713a123SBarry Smith 
19141713a123SBarry Smith 
19151713a123SBarry Smith 
1916