xref: /petsc/src/ts/interface/ts.c (revision d763cef207eb7fa01b459cfa1b5fffdd1cf73419)
1*d763cef2SBarry Smith #ifdef PETSC_RCS_HEADER
2*d763cef2SBarry Smith 
3*d763cef2SBarry Smith #endif
4*d763cef2SBarry Smith 
5*d763cef2SBarry Smith #include "src/ts/tsimpl.h"        /*I "ts.h"  I*/
6*d763cef2SBarry Smith 
7*d763cef2SBarry Smith #undef __FUNC__
8*d763cef2SBarry Smith #define __FUNC__ "TSComputeRHSFunction"
9*d763cef2SBarry Smith /*
10*d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
11*d763cef2SBarry Smith 
12*d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
13*d763cef2SBarry Smith    this routine applies the matrix.
14*d763cef2SBarry Smith */
15*d763cef2SBarry Smith int TSComputeRHSFunction(TS ts,double t,Vec x, Vec y)
16*d763cef2SBarry Smith {
17*d763cef2SBarry Smith   int ierr;
18*d763cef2SBarry Smith 
19*d763cef2SBarry Smith   PetscFunctionBegin;
20*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
21*d763cef2SBarry Smith   PetscValidHeader(x);  PetscValidHeader(y);
22*d763cef2SBarry Smith 
23*d763cef2SBarry Smith   if (ts->rhsfunction) {
24*d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
25*d763cef2SBarry Smith     ierr = (*ts->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
26*d763cef2SBarry Smith     PetscStackPop;
27*d763cef2SBarry Smith     PetscFunctionReturn(0);
28*d763cef2SBarry Smith   }
29*d763cef2SBarry Smith 
30*d763cef2SBarry Smith   if (ts->rhsmatrix) { /* assemble matrix for this timestep */
31*d763cef2SBarry Smith     MatStructure flg;
32*d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side matrix function");
33*d763cef2SBarry Smith     ierr = (*ts->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP); CHKERRQ(ierr);
34*d763cef2SBarry Smith     PetscStackPop;
35*d763cef2SBarry Smith   }
36*d763cef2SBarry Smith   ierr = MatMult(ts->A,x,y); CHKERRQ(ierr);
37*d763cef2SBarry Smith 
38*d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
39*d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y); CHKERRQ(ierr);
40*d763cef2SBarry Smith 
41*d763cef2SBarry Smith   PetscFunctionReturn(0);
42*d763cef2SBarry Smith }
43*d763cef2SBarry Smith 
44*d763cef2SBarry Smith #undef __FUNC__
45*d763cef2SBarry Smith #define __FUNC__ "TSSetRHSFunction"
46*d763cef2SBarry Smith /*@C
47*d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
48*d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
49*d763cef2SBarry Smith 
50*d763cef2SBarry Smith     Collective on TS
51*d763cef2SBarry Smith 
52*d763cef2SBarry Smith     Input Parameters:
53*d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
54*d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
55*d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
56*d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
57*d763cef2SBarry Smith 
58*d763cef2SBarry Smith     Calling sequence of func:
59*d763cef2SBarry Smith $     func (TS ts,double t,Vec u,Vec F,void *ctx);
60*d763cef2SBarry Smith 
61*d763cef2SBarry Smith +   t - current timestep
62*d763cef2SBarry Smith .   u - input vector
63*d763cef2SBarry Smith .   F - function vector
64*d763cef2SBarry Smith -   ctx - [optional] user-defined function context
65*d763cef2SBarry Smith 
66*d763cef2SBarry Smith     Important:
67*d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
68*d763cef2SBarry Smith 
69*d763cef2SBarry Smith     Level: beginner
70*d763cef2SBarry Smith 
71*d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
72*d763cef2SBarry Smith 
73*d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
74*d763cef2SBarry Smith @*/
75*d763cef2SBarry Smith int TSSetRHSFunction(TS ts,int (*f)(TS,double,Vec,Vec,void*),void *ctx)
76*d763cef2SBarry Smith {
77*d763cef2SBarry Smith   PetscFunctionBegin;
78*d763cef2SBarry Smith 
79*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
80*d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
81*d763cef2SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Cannot set function for linear problem");
82*d763cef2SBarry Smith   }
83*d763cef2SBarry Smith   ts->rhsfunction = f;
84*d763cef2SBarry Smith   ts->funP        = ctx;
85*d763cef2SBarry Smith   PetscFunctionReturn(0);
86*d763cef2SBarry Smith }
87*d763cef2SBarry Smith 
88*d763cef2SBarry Smith #undef __FUNC__
89*d763cef2SBarry Smith #define __FUNC__ "TSSetRHSMatrix"
90*d763cef2SBarry Smith /*@C
91*d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
92*d763cef2SBarry Smith    Also sets the location to store A.
93*d763cef2SBarry Smith 
94*d763cef2SBarry Smith    Collective on TS
95*d763cef2SBarry Smith 
96*d763cef2SBarry Smith    Input Parameters:
97*d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
98*d763cef2SBarry Smith .  A   - matrix
99*d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
100*d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
101*d763cef2SBarry Smith          if A is not a function of t.
102*d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
103*d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
104*d763cef2SBarry Smith 
105*d763cef2SBarry Smith    Calling sequence of func:
106*d763cef2SBarry Smith $     func (TS ts,double t,Mat *A,Mat *B,int *flag,void *ctx);
107*d763cef2SBarry Smith 
108*d763cef2SBarry Smith +  t - current timestep
109*d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
110*d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
111*d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
112*d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
113*d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
114*d763cef2SBarry Smith 
115*d763cef2SBarry Smith    Notes:
116*d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
117*d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
118*d763cef2SBarry Smith 
119*d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
120*d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
121*d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
122*d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
123*d763cef2SBarry Smith    throughout the global iterations.
124*d763cef2SBarry Smith 
125*d763cef2SBarry Smith    Important:
126*d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
127*d763cef2SBarry Smith 
128*d763cef2SBarry Smith    Level: beginner
129*d763cef2SBarry Smith 
130*d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
131*d763cef2SBarry Smith 
132*d763cef2SBarry Smith .seealso: TSSetRHSFunction()
133*d763cef2SBarry Smith @*/
134*d763cef2SBarry Smith int TSSetRHSMatrix(TS ts,Mat A, Mat B,int (*f)(TS,double,Mat*,Mat*,MatStructure*,void*),void *ctx)
135*d763cef2SBarry Smith {
136*d763cef2SBarry Smith   PetscFunctionBegin;
137*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
138*d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
139*d763cef2SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for nonlinear problems; use TSSetRHSJacobian()");
140*d763cef2SBarry Smith   }
141*d763cef2SBarry Smith 
142*d763cef2SBarry Smith   ts->rhsmatrix = f;
143*d763cef2SBarry Smith   ts->jacP      = ctx;
144*d763cef2SBarry Smith   ts->A         = A;
145*d763cef2SBarry Smith   ts->B         = B;
146*d763cef2SBarry Smith 
147*d763cef2SBarry Smith   PetscFunctionReturn(0);
148*d763cef2SBarry Smith }
149*d763cef2SBarry Smith 
150*d763cef2SBarry Smith #undef __FUNC__
151*d763cef2SBarry Smith #define __FUNC__ "TSSetRHSJacobian"
152*d763cef2SBarry Smith /*@C
153*d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
154*d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
155*d763cef2SBarry Smith 
156*d763cef2SBarry Smith    Collective on TS
157*d763cef2SBarry Smith 
158*d763cef2SBarry Smith    Input Parameters:
159*d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
160*d763cef2SBarry Smith .  A   - Jacobian matrix
161*d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
162*d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
163*d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
164*d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
165*d763cef2SBarry Smith 
166*d763cef2SBarry Smith    Calling sequence of func:
167*d763cef2SBarry Smith $     func (TS ts,double t,Vec u,Mat *A,Mat *B,int *flag,void *ctx);
168*d763cef2SBarry Smith 
169*d763cef2SBarry Smith +  t - current timestep
170*d763cef2SBarry Smith .  u - input vector
171*d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
172*d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
173*d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
174*d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
175*d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
176*d763cef2SBarry Smith 
177*d763cef2SBarry Smith    Notes:
178*d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
179*d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
180*d763cef2SBarry Smith 
181*d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
182*d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
183*d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
184*d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
185*d763cef2SBarry Smith    throughout the global iterations.
186*d763cef2SBarry Smith 
187*d763cef2SBarry Smith    Level: beginner
188*d763cef2SBarry Smith 
189*d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
190*d763cef2SBarry Smith 
191*d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
192*d763cef2SBarry Smith           SNESDefaultComputeJacobianColor()
193*d763cef2SBarry Smith 
194*d763cef2SBarry Smith @*/
195*d763cef2SBarry Smith int TSSetRHSJacobian(TS ts,Mat A, Mat B,int (*f)(TS,double,Vec,Mat*,Mat*,
196*d763cef2SBarry Smith                      MatStructure*,void*),void *ctx)
197*d763cef2SBarry Smith {
198*d763cef2SBarry Smith   PetscFunctionBegin;
199*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
200*d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
201*d763cef2SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for linear problems; use TSSetRHSMatrix()");
202*d763cef2SBarry Smith   }
203*d763cef2SBarry Smith 
204*d763cef2SBarry Smith   ts->rhsjacobian = f;
205*d763cef2SBarry Smith   ts->jacP        = ctx;
206*d763cef2SBarry Smith   ts->A           = A;
207*d763cef2SBarry Smith   ts->B           = B;
208*d763cef2SBarry Smith   PetscFunctionReturn(0);
209*d763cef2SBarry Smith }
210*d763cef2SBarry Smith 
211*d763cef2SBarry Smith #undef __FUNC__
212*d763cef2SBarry Smith #define __FUNC__ "TSComputeRHSBoundaryConditions"
213*d763cef2SBarry Smith /*
214*d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
215*d763cef2SBarry Smith 
216*d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
217*d763cef2SBarry Smith    this routine applies the matrix.
218*d763cef2SBarry Smith */
219*d763cef2SBarry Smith int TSComputeRHSBoundaryConditions(TS ts,double t,Vec x)
220*d763cef2SBarry Smith {
221*d763cef2SBarry Smith   int ierr;
222*d763cef2SBarry Smith 
223*d763cef2SBarry Smith   PetscFunctionBegin;
224*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
225*d763cef2SBarry Smith   PetscValidHeader(x);
226*d763cef2SBarry Smith 
227*d763cef2SBarry Smith   if (ts->rhsbc) {
228*d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
229*d763cef2SBarry Smith     ierr = (*ts->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
230*d763cef2SBarry Smith     PetscStackPop;
231*d763cef2SBarry Smith     PetscFunctionReturn(0);
232*d763cef2SBarry Smith   }
233*d763cef2SBarry Smith 
234*d763cef2SBarry Smith   PetscFunctionReturn(0);
235*d763cef2SBarry Smith }
236*d763cef2SBarry Smith 
237*d763cef2SBarry Smith #undef __FUNC__
238*d763cef2SBarry Smith #define __FUNC__ "TSSetRHSBoundaryConditions"
239*d763cef2SBarry Smith /*@C
240*d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
241*d763cef2SBarry Smith     boundary conditions for the function F.
242*d763cef2SBarry Smith 
243*d763cef2SBarry Smith     Collective on TS
244*d763cef2SBarry Smith 
245*d763cef2SBarry Smith     Input Parameters:
246*d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
247*d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
248*d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
249*d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
250*d763cef2SBarry Smith 
251*d763cef2SBarry Smith     Calling sequence of func:
252*d763cef2SBarry Smith $     func (TS ts,double t,Vec F,void *ctx);
253*d763cef2SBarry Smith 
254*d763cef2SBarry Smith +   t - current timestep
255*d763cef2SBarry Smith .   F - function vector
256*d763cef2SBarry Smith -   ctx - [optional] user-defined function context
257*d763cef2SBarry Smith 
258*d763cef2SBarry Smith     Level: intermediate
259*d763cef2SBarry Smith 
260*d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
261*d763cef2SBarry Smith @*/
262*d763cef2SBarry Smith int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,double,Vec,void*),void *ctx)
263*d763cef2SBarry Smith {
264*d763cef2SBarry Smith   PetscFunctionBegin;
265*d763cef2SBarry Smith 
266*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
267*d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
268*d763cef2SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For linear problems only");
269*d763cef2SBarry Smith   }
270*d763cef2SBarry Smith   ts->rhsbc = f;
271*d763cef2SBarry Smith   ts->bcP   = ctx;
272*d763cef2SBarry Smith   PetscFunctionReturn(0);
273*d763cef2SBarry Smith }
274*d763cef2SBarry Smith 
275*d763cef2SBarry Smith #undef __FUNC__
276*d763cef2SBarry Smith #define __FUNC__ "TSView"
277*d763cef2SBarry Smith /*@
278*d763cef2SBarry Smith     TSView - Prints the TS data structure.
279*d763cef2SBarry Smith 
280*d763cef2SBarry Smith     Collective on TS, unless Viewer is VIEWER_STDOUT_SELF
281*d763cef2SBarry Smith 
282*d763cef2SBarry Smith     Input Parameters:
283*d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
284*d763cef2SBarry Smith -   viewer - visualization context
285*d763cef2SBarry Smith 
286*d763cef2SBarry Smith     Options Database Key:
287*d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
288*d763cef2SBarry Smith 
289*d763cef2SBarry Smith     Notes:
290*d763cef2SBarry Smith     The available visualization contexts include
291*d763cef2SBarry Smith +     VIEWER_STDOUT_SELF - standard output (default)
292*d763cef2SBarry Smith -     VIEWER_STDOUT_WORLD - synchronized standard
293*d763cef2SBarry Smith          output where only the first processor opens
294*d763cef2SBarry Smith          the file.  All other processors send their
295*d763cef2SBarry Smith          data to the first processor to print.
296*d763cef2SBarry Smith 
297*d763cef2SBarry Smith     The user can open an alternative visualization context with
298*d763cef2SBarry Smith     ViewerASCIIOpen() - output to a specified file.
299*d763cef2SBarry Smith 
300*d763cef2SBarry Smith     Level: beginner
301*d763cef2SBarry Smith 
302*d763cef2SBarry Smith .keywords: TS, timestep, view
303*d763cef2SBarry Smith 
304*d763cef2SBarry Smith .seealso: ViewerASCIIOpen()
305*d763cef2SBarry Smith @*/
306*d763cef2SBarry Smith int TSView(TS ts,Viewer viewer)
307*d763cef2SBarry Smith {
308*d763cef2SBarry Smith   int                 ierr;
309*d763cef2SBarry Smith   char                *method;
310*d763cef2SBarry Smith   ViewerType          vtype;
311*d763cef2SBarry Smith 
312*d763cef2SBarry Smith   PetscFunctionBegin;
313*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
314*d763cef2SBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
315*d763cef2SBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
316*d763cef2SBarry Smith     ViewerASCIIPrintf(viewer,"TS Object:\n");
317*d763cef2SBarry Smith     ierr = TSGetType(ts,(TSType *)&method);CHKERRQ(ierr);
318*d763cef2SBarry Smith     ViewerASCIIPrintf(viewer,"  method: %s\n",method);
319*d763cef2SBarry Smith     if (ts->view) {
320*d763cef2SBarry Smith       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
321*d763cef2SBarry Smith       ierr = (*ts->view)(ts,viewer);CHKERRQ(ierr);
322*d763cef2SBarry Smith       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
323*d763cef2SBarry Smith     }
324*d763cef2SBarry Smith     ViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);
325*d763cef2SBarry Smith     ViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);
326*d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
327*d763cef2SBarry Smith       ViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);
328*d763cef2SBarry Smith     }
329*d763cef2SBarry Smith     ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);
330*d763cef2SBarry Smith   } else if (PetscTypeCompare(vtype,STRING_VIEWER)) {
331*d763cef2SBarry Smith     ierr = TSGetType(ts,(TSType *)&method);CHKERRQ(ierr);
332*d763cef2SBarry Smith     ierr = ViewerStringSPrintf(viewer," %-7.7s",method);CHKERRQ(ierr);
333*d763cef2SBarry Smith   }
334*d763cef2SBarry Smith   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
335*d763cef2SBarry Smith   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
336*d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
337*d763cef2SBarry Smith   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
338*d763cef2SBarry Smith   PetscFunctionReturn(0);
339*d763cef2SBarry Smith }
340*d763cef2SBarry Smith 
341*d763cef2SBarry Smith 
342*d763cef2SBarry Smith #undef __FUNC__
343*d763cef2SBarry Smith #define __FUNC__ "TSSetApplicationContext"
344*d763cef2SBarry Smith /*@C
345*d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
346*d763cef2SBarry Smith    the timesteppers.
347*d763cef2SBarry Smith 
348*d763cef2SBarry Smith    Collective on TS
349*d763cef2SBarry Smith 
350*d763cef2SBarry Smith    Input Parameters:
351*d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
352*d763cef2SBarry Smith -  usrP - optional user context
353*d763cef2SBarry Smith 
354*d763cef2SBarry Smith    Level: intermediate
355*d763cef2SBarry Smith 
356*d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
357*d763cef2SBarry Smith 
358*d763cef2SBarry Smith .seealso: TSGetApplicationContext()
359*d763cef2SBarry Smith @*/
360*d763cef2SBarry Smith int TSSetApplicationContext(TS ts,void *usrP)
361*d763cef2SBarry Smith {
362*d763cef2SBarry Smith   PetscFunctionBegin;
363*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
364*d763cef2SBarry Smith   ts->user = usrP;
365*d763cef2SBarry Smith   PetscFunctionReturn(0);
366*d763cef2SBarry Smith }
367*d763cef2SBarry Smith 
368*d763cef2SBarry Smith #undef __FUNC__
369*d763cef2SBarry Smith #define __FUNC__ "TSGetApplicationContext"
370*d763cef2SBarry Smith /*@C
371*d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
372*d763cef2SBarry Smith     timestepper.
373*d763cef2SBarry Smith 
374*d763cef2SBarry Smith     Not Collective
375*d763cef2SBarry Smith 
376*d763cef2SBarry Smith     Input Parameter:
377*d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
378*d763cef2SBarry Smith 
379*d763cef2SBarry Smith     Output Parameter:
380*d763cef2SBarry Smith .   usrP - user context
381*d763cef2SBarry Smith 
382*d763cef2SBarry Smith     Level: intermediate
383*d763cef2SBarry Smith 
384*d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
385*d763cef2SBarry Smith 
386*d763cef2SBarry Smith .seealso: TSSetApplicationContext()
387*d763cef2SBarry Smith @*/
388*d763cef2SBarry Smith int TSGetApplicationContext( TS ts,  void **usrP )
389*d763cef2SBarry Smith {
390*d763cef2SBarry Smith   PetscFunctionBegin;
391*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
392*d763cef2SBarry Smith   *usrP = ts->user;
393*d763cef2SBarry Smith   PetscFunctionReturn(0);
394*d763cef2SBarry Smith }
395*d763cef2SBarry Smith 
396*d763cef2SBarry Smith #undef __FUNC__
397*d763cef2SBarry Smith #define __FUNC__ "TSGetTimeStepNumber"
398*d763cef2SBarry Smith /*@
399*d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
400*d763cef2SBarry Smith 
401*d763cef2SBarry Smith    Not Collective
402*d763cef2SBarry Smith 
403*d763cef2SBarry Smith    Input Parameter:
404*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
405*d763cef2SBarry Smith 
406*d763cef2SBarry Smith    Output Parameter:
407*d763cef2SBarry Smith .  iter - number steps so far
408*d763cef2SBarry Smith 
409*d763cef2SBarry Smith    Level: intermediate
410*d763cef2SBarry Smith 
411*d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
412*d763cef2SBarry Smith @*/
413*d763cef2SBarry Smith int TSGetTimeStepNumber(TS ts,int* iter)
414*d763cef2SBarry Smith {
415*d763cef2SBarry Smith   PetscFunctionBegin;
416*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
417*d763cef2SBarry Smith   *iter = ts->steps;
418*d763cef2SBarry Smith   PetscFunctionReturn(0);
419*d763cef2SBarry Smith }
420*d763cef2SBarry Smith 
421*d763cef2SBarry Smith #undef __FUNC__
422*d763cef2SBarry Smith #define __FUNC__ "TSSetInitialTimeStep"
423*d763cef2SBarry Smith /*@
424*d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
425*d763cef2SBarry Smith    as well as the initial time.
426*d763cef2SBarry Smith 
427*d763cef2SBarry Smith    Collective on TS
428*d763cef2SBarry Smith 
429*d763cef2SBarry Smith    Input Parameters:
430*d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
431*d763cef2SBarry Smith .  initial_time - the initial time
432*d763cef2SBarry Smith -  time_step - the size of the timestep
433*d763cef2SBarry Smith 
434*d763cef2SBarry Smith    Level: intermediate
435*d763cef2SBarry Smith 
436*d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
437*d763cef2SBarry Smith 
438*d763cef2SBarry Smith .keywords: TS, set, initial, timestep
439*d763cef2SBarry Smith @*/
440*d763cef2SBarry Smith int TSSetInitialTimeStep(TS ts,double initial_time,double time_step)
441*d763cef2SBarry Smith {
442*d763cef2SBarry Smith   PetscFunctionBegin;
443*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
444*d763cef2SBarry Smith   ts->time_step         = time_step;
445*d763cef2SBarry Smith   ts->initial_time_step = time_step;
446*d763cef2SBarry Smith   ts->ptime             = initial_time;
447*d763cef2SBarry Smith   PetscFunctionReturn(0);
448*d763cef2SBarry Smith }
449*d763cef2SBarry Smith 
450*d763cef2SBarry Smith #undef __FUNC__
451*d763cef2SBarry Smith #define __FUNC__ "TSSetTimeStep"
452*d763cef2SBarry Smith /*@
453*d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
454*d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
455*d763cef2SBarry Smith 
456*d763cef2SBarry Smith    Collective on TS
457*d763cef2SBarry Smith 
458*d763cef2SBarry Smith    Input Parameters:
459*d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
460*d763cef2SBarry Smith -  time_step - the size of the timestep
461*d763cef2SBarry Smith 
462*d763cef2SBarry Smith    Level: intermediate
463*d763cef2SBarry Smith 
464*d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
465*d763cef2SBarry Smith 
466*d763cef2SBarry Smith .keywords: TS, set, timestep
467*d763cef2SBarry Smith @*/
468*d763cef2SBarry Smith int TSSetTimeStep(TS ts,double time_step)
469*d763cef2SBarry Smith {
470*d763cef2SBarry Smith   PetscFunctionBegin;
471*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
472*d763cef2SBarry Smith   ts->time_step = time_step;
473*d763cef2SBarry Smith   PetscFunctionReturn(0);
474*d763cef2SBarry Smith }
475*d763cef2SBarry Smith 
476*d763cef2SBarry Smith #undef __FUNC__
477*d763cef2SBarry Smith #define __FUNC__ "TSGetTimeStep"
478*d763cef2SBarry Smith /*@
479*d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
480*d763cef2SBarry Smith 
481*d763cef2SBarry Smith    Not Collective
482*d763cef2SBarry Smith 
483*d763cef2SBarry Smith    Input Parameter:
484*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
485*d763cef2SBarry Smith 
486*d763cef2SBarry Smith    Output Parameter:
487*d763cef2SBarry Smith .  dt - the current timestep size
488*d763cef2SBarry Smith 
489*d763cef2SBarry Smith    Level: intermediate
490*d763cef2SBarry Smith 
491*d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
492*d763cef2SBarry Smith 
493*d763cef2SBarry Smith .keywords: TS, get, timestep
494*d763cef2SBarry Smith @*/
495*d763cef2SBarry Smith int TSGetTimeStep(TS ts,double* dt)
496*d763cef2SBarry Smith {
497*d763cef2SBarry Smith   PetscFunctionBegin;
498*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
499*d763cef2SBarry Smith   *dt = ts->time_step;
500*d763cef2SBarry Smith   PetscFunctionReturn(0);
501*d763cef2SBarry Smith }
502*d763cef2SBarry Smith 
503*d763cef2SBarry Smith #undef __FUNC__
504*d763cef2SBarry Smith #define __FUNC__ "TSGetSolution"
505*d763cef2SBarry Smith /*@C
506*d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
507*d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
508*d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
509*d763cef2SBarry Smith    the solution at the next timestep has been calculated.
510*d763cef2SBarry Smith 
511*d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
512*d763cef2SBarry Smith 
513*d763cef2SBarry Smith    Input Parameter:
514*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
515*d763cef2SBarry Smith 
516*d763cef2SBarry Smith    Output Parameter:
517*d763cef2SBarry Smith .  v - the vector containing the solution
518*d763cef2SBarry Smith 
519*d763cef2SBarry Smith    Level: intermediate
520*d763cef2SBarry Smith 
521*d763cef2SBarry Smith .seealso: TSGetTimeStep()
522*d763cef2SBarry Smith 
523*d763cef2SBarry Smith .keywords: TS, timestep, get, solution
524*d763cef2SBarry Smith @*/
525*d763cef2SBarry Smith int TSGetSolution(TS ts,Vec *v)
526*d763cef2SBarry Smith {
527*d763cef2SBarry Smith   PetscFunctionBegin;
528*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
529*d763cef2SBarry Smith   *v = ts->vec_sol_always;
530*d763cef2SBarry Smith   PetscFunctionReturn(0);
531*d763cef2SBarry Smith }
532*d763cef2SBarry Smith 
533*d763cef2SBarry Smith #undef __FUNC__
534*d763cef2SBarry Smith #define __FUNC__ "TSPublish_Petsc"
535*d763cef2SBarry Smith static int TSPublish_Petsc(PetscObject object)
536*d763cef2SBarry Smith {
537*d763cef2SBarry Smith #if defined(HAVE_AMS)
538*d763cef2SBarry Smith   TS   v = (TS) object;
539*d763cef2SBarry Smith   int  ierr;
540*d763cef2SBarry Smith 
541*d763cef2SBarry Smith   PetscFunctionBegin;
542*d763cef2SBarry Smith 
543*d763cef2SBarry Smith   /* if it is already published then return */
544*d763cef2SBarry Smith   if (v->amem >=0 ) PetscFunctionReturn(0);
545*d763cef2SBarry Smith 
546*d763cef2SBarry Smith   ierr = PetscObjectPublishBaseBegin(object);CHKERRQ(ierr);
547*d763cef2SBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ,
548*d763cef2SBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
549*d763cef2SBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ,
550*d763cef2SBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
551*d763cef2SBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1,
552*d763cef2SBarry Smith                                AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
553*d763cef2SBarry Smith   ierr = PetscObjectPublishBaseEnd(object);CHKERRQ(ierr);
554*d763cef2SBarry Smith #else
555*d763cef2SBarry Smith   PetscFunctionBegin;
556*d763cef2SBarry Smith #endif
557*d763cef2SBarry Smith   PetscFunctionReturn(0);
558*d763cef2SBarry Smith }
559*d763cef2SBarry Smith 
560*d763cef2SBarry Smith /* -----------------------------------------------------------*/
561*d763cef2SBarry Smith 
562*d763cef2SBarry Smith #undef __FUNC__
563*d763cef2SBarry Smith #define __FUNC__ "TSCreate"
564*d763cef2SBarry Smith /*@C
565*d763cef2SBarry Smith    TSCreate - Creates a timestepper context.
566*d763cef2SBarry Smith 
567*d763cef2SBarry Smith    Collective on MPI_Comm
568*d763cef2SBarry Smith 
569*d763cef2SBarry Smith    Input Parameter:
570*d763cef2SBarry Smith +  comm - MPI communicator
571*d763cef2SBarry Smith -  type - One of  TS_LINEAR,TS_NONLINEAR
572*d763cef2SBarry Smith    where these types refer to problems of the forms
573*d763cef2SBarry Smith .vb
574*d763cef2SBarry Smith          U_t = A U
575*d763cef2SBarry Smith          U_t = A(t) U
576*d763cef2SBarry Smith          U_t = F(t,U)
577*d763cef2SBarry Smith .ve
578*d763cef2SBarry Smith 
579*d763cef2SBarry Smith    Output Parameter:
580*d763cef2SBarry Smith .  outts - the new TS context
581*d763cef2SBarry Smith 
582*d763cef2SBarry Smith    Level: beginner
583*d763cef2SBarry Smith 
584*d763cef2SBarry Smith .keywords: TS, timestep, create, context
585*d763cef2SBarry Smith 
586*d763cef2SBarry Smith .seealso: TSSetUp(), TSStep(), TSDestroy()
587*d763cef2SBarry Smith @*/
588*d763cef2SBarry Smith int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts)
589*d763cef2SBarry Smith {
590*d763cef2SBarry Smith   TS   ts;
591*d763cef2SBarry Smith 
592*d763cef2SBarry Smith   PetscFunctionBegin;
593*d763cef2SBarry Smith   *outts = 0;
594*d763cef2SBarry Smith   PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView);
595*d763cef2SBarry Smith   PLogObjectCreate(ts);
596*d763cef2SBarry Smith   ts->bops->publish     = TSPublish_Petsc;
597*d763cef2SBarry Smith   ts->max_steps         = 5000;
598*d763cef2SBarry Smith   ts->max_time          = 5.0;
599*d763cef2SBarry Smith   ts->time_step         = .1;
600*d763cef2SBarry Smith   ts->initial_time_step = ts->time_step;
601*d763cef2SBarry Smith   ts->steps             = 0;
602*d763cef2SBarry Smith   ts->ptime             = 0.0;
603*d763cef2SBarry Smith   ts->data              = 0;
604*d763cef2SBarry Smith   ts->view              = 0;
605*d763cef2SBarry Smith   ts->setupcalled       = 0;
606*d763cef2SBarry Smith   ts->problem_type      = problemtype;
607*d763cef2SBarry Smith   ts->numbermonitors    = 0;
608*d763cef2SBarry Smith   ts->linear_its        = 0;
609*d763cef2SBarry Smith   ts->nonlinear_its     = 0;
610*d763cef2SBarry Smith 
611*d763cef2SBarry Smith   *outts = ts;
612*d763cef2SBarry Smith   PetscFunctionReturn(0);
613*d763cef2SBarry Smith }
614*d763cef2SBarry Smith 
615*d763cef2SBarry Smith /* ----- Routines to initialize and destroy a timestepper ---- */
616*d763cef2SBarry Smith 
617*d763cef2SBarry Smith #undef __FUNC__
618*d763cef2SBarry Smith #define __FUNC__ "TSSetUp"
619*d763cef2SBarry Smith /*@
620*d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
621*d763cef2SBarry Smith    of a timestepper.
622*d763cef2SBarry Smith 
623*d763cef2SBarry Smith    Collective on TS
624*d763cef2SBarry Smith 
625*d763cef2SBarry Smith    Input Parameter:
626*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
627*d763cef2SBarry Smith 
628*d763cef2SBarry Smith    Notes:
629*d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
630*d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
631*d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
632*d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
633*d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
634*d763cef2SBarry Smith 
635*d763cef2SBarry Smith    Level: advanced
636*d763cef2SBarry Smith 
637*d763cef2SBarry Smith .keywords: TS, timestep, setup
638*d763cef2SBarry Smith 
639*d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
640*d763cef2SBarry Smith @*/
641*d763cef2SBarry Smith int TSSetUp(TS ts)
642*d763cef2SBarry Smith {
643*d763cef2SBarry Smith   int ierr;
644*d763cef2SBarry Smith 
645*d763cef2SBarry Smith   PetscFunctionBegin;
646*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
647*d763cef2SBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call TSSetSolution() first");
648*d763cef2SBarry Smith   if (!ts->type_name) {
649*d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
650*d763cef2SBarry Smith   }
651*d763cef2SBarry Smith   ierr = (*ts->setup)(ts); CHKERRQ(ierr);
652*d763cef2SBarry Smith   ts->setupcalled = 1;
653*d763cef2SBarry Smith   PetscFunctionReturn(0);
654*d763cef2SBarry Smith }
655*d763cef2SBarry Smith 
656*d763cef2SBarry Smith #undef __FUNC__
657*d763cef2SBarry Smith #define __FUNC__ "TSDestroy"
658*d763cef2SBarry Smith /*@C
659*d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
660*d763cef2SBarry Smith    with TSCreate().
661*d763cef2SBarry Smith 
662*d763cef2SBarry Smith    Collective on TS
663*d763cef2SBarry Smith 
664*d763cef2SBarry Smith    Input Parameter:
665*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
666*d763cef2SBarry Smith 
667*d763cef2SBarry Smith    Level: beginner
668*d763cef2SBarry Smith 
669*d763cef2SBarry Smith .keywords: TS, timestepper, destroy
670*d763cef2SBarry Smith 
671*d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
672*d763cef2SBarry Smith @*/
673*d763cef2SBarry Smith int TSDestroy(TS ts)
674*d763cef2SBarry Smith {
675*d763cef2SBarry Smith   int ierr;
676*d763cef2SBarry Smith 
677*d763cef2SBarry Smith   PetscFunctionBegin;
678*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
679*d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
680*d763cef2SBarry Smith 
681*d763cef2SBarry Smith   if (ts->sles) {ierr = SLESDestroy(ts->sles); CHKERRQ(ierr);}
682*d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes); CHKERRQ(ierr);}
683*d763cef2SBarry Smith   ierr = (*(ts)->destroy)(ts); CHKERRQ(ierr);
684*d763cef2SBarry Smith   PLogObjectDestroy((PetscObject)ts);
685*d763cef2SBarry Smith   PetscHeaderDestroy((PetscObject)ts);
686*d763cef2SBarry Smith   PetscFunctionReturn(0);
687*d763cef2SBarry Smith }
688*d763cef2SBarry Smith 
689*d763cef2SBarry Smith #undef __FUNC__
690*d763cef2SBarry Smith #define __FUNC__ "TSGetSNES"
691*d763cef2SBarry Smith /*@C
692*d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
693*d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
694*d763cef2SBarry Smith 
695*d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
696*d763cef2SBarry Smith 
697*d763cef2SBarry Smith    Input Parameter:
698*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
699*d763cef2SBarry Smith 
700*d763cef2SBarry Smith    Output Parameter:
701*d763cef2SBarry Smith .  snes - the nonlinear solver context
702*d763cef2SBarry Smith 
703*d763cef2SBarry Smith    Notes:
704*d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
705*d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
706*d763cef2SBarry Smith    SLES, KSP, and PC contexts as well.
707*d763cef2SBarry Smith 
708*d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
709*d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
710*d763cef2SBarry Smith 
711*d763cef2SBarry Smith    Level: beginner
712*d763cef2SBarry Smith 
713*d763cef2SBarry Smith .keywords: timestep, get, SNES
714*d763cef2SBarry Smith @*/
715*d763cef2SBarry Smith int TSGetSNES(TS ts,SNES *snes)
716*d763cef2SBarry Smith {
717*d763cef2SBarry Smith   PetscFunctionBegin;
718*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
719*d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Nonlinear only; use TSGetSLES()");
720*d763cef2SBarry Smith   *snes = ts->snes;
721*d763cef2SBarry Smith   PetscFunctionReturn(0);
722*d763cef2SBarry Smith }
723*d763cef2SBarry Smith 
724*d763cef2SBarry Smith #undef __FUNC__
725*d763cef2SBarry Smith #define __FUNC__ "TSGetSLES"
726*d763cef2SBarry Smith /*@C
727*d763cef2SBarry Smith    TSGetSLES - Returns the SLES (linear solver) associated with
728*d763cef2SBarry Smith    a TS (timestepper) context.
729*d763cef2SBarry Smith 
730*d763cef2SBarry Smith    Not Collective, but SLES is parallel if TS is parallel
731*d763cef2SBarry Smith 
732*d763cef2SBarry Smith    Input Parameter:
733*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
734*d763cef2SBarry Smith 
735*d763cef2SBarry Smith    Output Parameter:
736*d763cef2SBarry Smith .  sles - the nonlinear solver context
737*d763cef2SBarry Smith 
738*d763cef2SBarry Smith    Notes:
739*d763cef2SBarry Smith    The user can then directly manipulate the SLES context to set various
740*d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
741*d763cef2SBarry Smith    KSP and PC contexts as well.
742*d763cef2SBarry Smith 
743*d763cef2SBarry Smith    TSGetSLES() does not work for integrators that do not use SLES;
744*d763cef2SBarry Smith    in this case TSGetSLES() returns PETSC_NULL in sles.
745*d763cef2SBarry Smith 
746*d763cef2SBarry Smith    Level: beginner
747*d763cef2SBarry Smith 
748*d763cef2SBarry Smith .keywords: timestep, get, SLES
749*d763cef2SBarry Smith @*/
750*d763cef2SBarry Smith int TSGetSLES(TS ts,SLES *sles)
751*d763cef2SBarry Smith {
752*d763cef2SBarry Smith   PetscFunctionBegin;
753*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
754*d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Linear only; use TSGetSNES()");
755*d763cef2SBarry Smith   *sles = ts->sles;
756*d763cef2SBarry Smith   PetscFunctionReturn(0);
757*d763cef2SBarry Smith }
758*d763cef2SBarry Smith 
759*d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
760*d763cef2SBarry Smith 
761*d763cef2SBarry Smith #undef __FUNC__
762*d763cef2SBarry Smith #define __FUNC__ "TSSetDuration"
763*d763cef2SBarry Smith /*@
764*d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
765*d763cef2SBarry Smith    maximum time for iteration.
766*d763cef2SBarry Smith 
767*d763cef2SBarry Smith    Collective on TS
768*d763cef2SBarry Smith 
769*d763cef2SBarry Smith    Input Parameters:
770*d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
771*d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
772*d763cef2SBarry Smith -  maxtime - final time to iterate to
773*d763cef2SBarry Smith 
774*d763cef2SBarry Smith    Options Database Keys:
775*d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
776*d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
777*d763cef2SBarry Smith 
778*d763cef2SBarry Smith    Notes:
779*d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
780*d763cef2SBarry Smith 
781*d763cef2SBarry Smith    Level: intermediate
782*d763cef2SBarry Smith 
783*d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
784*d763cef2SBarry Smith @*/
785*d763cef2SBarry Smith int TSSetDuration(TS ts,int maxsteps,double maxtime)
786*d763cef2SBarry Smith {
787*d763cef2SBarry Smith   PetscFunctionBegin;
788*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
789*d763cef2SBarry Smith   ts->max_steps = maxsteps;
790*d763cef2SBarry Smith   ts->max_time  = maxtime;
791*d763cef2SBarry Smith   PetscFunctionReturn(0);
792*d763cef2SBarry Smith }
793*d763cef2SBarry Smith 
794*d763cef2SBarry Smith #undef __FUNC__
795*d763cef2SBarry Smith #define __FUNC__ "TSSetSolution"
796*d763cef2SBarry Smith /*@
797*d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
798*d763cef2SBarry Smith    for use by the TS routines.
799*d763cef2SBarry Smith 
800*d763cef2SBarry Smith    Collective on TS and Vec
801*d763cef2SBarry Smith 
802*d763cef2SBarry Smith    Input Parameters:
803*d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
804*d763cef2SBarry Smith -  x - the solution vector
805*d763cef2SBarry Smith 
806*d763cef2SBarry Smith    Level: beginner
807*d763cef2SBarry Smith 
808*d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
809*d763cef2SBarry Smith @*/
810*d763cef2SBarry Smith int TSSetSolution(TS ts,Vec x)
811*d763cef2SBarry Smith {
812*d763cef2SBarry Smith   PetscFunctionBegin;
813*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
814*d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
815*d763cef2SBarry Smith   PetscFunctionReturn(0);
816*d763cef2SBarry Smith }
817*d763cef2SBarry Smith 
818*d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
819*d763cef2SBarry Smith 
820*d763cef2SBarry Smith #undef __FUNC__
821*d763cef2SBarry Smith #define __FUNC__ "TSSetMonitor"
822*d763cef2SBarry Smith /*@C
823*d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
824*d763cef2SBarry Smith    timestep to display the iteration's  progress.
825*d763cef2SBarry Smith 
826*d763cef2SBarry Smith    Collective on TS
827*d763cef2SBarry Smith 
828*d763cef2SBarry Smith    Input Parameters:
829*d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
830*d763cef2SBarry Smith .  func - monitoring routine
831*d763cef2SBarry Smith -  mctx - [optional] user-defined context for private data for the
832*d763cef2SBarry Smith           monitor routine (may be PETSC_NULL)
833*d763cef2SBarry Smith 
834*d763cef2SBarry Smith    Calling sequence of func:
835*d763cef2SBarry Smith $    int func(TS ts,int steps,double time,Vec x,void *mctx)
836*d763cef2SBarry Smith 
837*d763cef2SBarry Smith +    ts - the TS context
838*d763cef2SBarry Smith .    steps - iteration number
839*d763cef2SBarry Smith .    time - current timestep
840*d763cef2SBarry Smith .    x - current iterate
841*d763cef2SBarry Smith -    mctx - [optional] monitoring context
842*d763cef2SBarry Smith 
843*d763cef2SBarry Smith    Notes:
844*d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
845*d763cef2SBarry Smith    already has been loaded.
846*d763cef2SBarry Smith 
847*d763cef2SBarry Smith    Level: intermediate
848*d763cef2SBarry Smith 
849*d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
850*d763cef2SBarry Smith 
851*d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
852*d763cef2SBarry Smith @*/
853*d763cef2SBarry Smith int TSSetMonitor(TS ts, int (*monitor)(TS,int,double,Vec,void*), void *mctx )
854*d763cef2SBarry Smith {
855*d763cef2SBarry Smith   PetscFunctionBegin;
856*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
857*d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
858*d763cef2SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
859*d763cef2SBarry Smith   }
860*d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
861*d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
862*d763cef2SBarry Smith   PetscFunctionReturn(0);
863*d763cef2SBarry Smith }
864*d763cef2SBarry Smith 
865*d763cef2SBarry Smith #undef __FUNC__
866*d763cef2SBarry Smith #define __FUNC__ "TSClearMonitor"
867*d763cef2SBarry Smith /*@C
868*d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
869*d763cef2SBarry Smith 
870*d763cef2SBarry Smith    Collective on TS
871*d763cef2SBarry Smith 
872*d763cef2SBarry Smith    Input Parameters:
873*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
874*d763cef2SBarry Smith 
875*d763cef2SBarry Smith    Notes:
876*d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
877*d763cef2SBarry Smith 
878*d763cef2SBarry Smith    Level: intermediate
879*d763cef2SBarry Smith 
880*d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
881*d763cef2SBarry Smith 
882*d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
883*d763cef2SBarry Smith @*/
884*d763cef2SBarry Smith int TSClearMonitor(TS ts)
885*d763cef2SBarry Smith {
886*d763cef2SBarry Smith   PetscFunctionBegin;
887*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
888*d763cef2SBarry Smith   ts->numbermonitors = 0;
889*d763cef2SBarry Smith   PetscFunctionReturn(0);
890*d763cef2SBarry Smith }
891*d763cef2SBarry Smith 
892*d763cef2SBarry Smith #undef __FUNC__
893*d763cef2SBarry Smith #define __FUNC__ "TSDefaultMonitor"
894*d763cef2SBarry Smith int TSDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx)
895*d763cef2SBarry Smith {
896*d763cef2SBarry Smith   PetscFunctionBegin;
897*d763cef2SBarry Smith   PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,time);
898*d763cef2SBarry Smith   PetscFunctionReturn(0);
899*d763cef2SBarry Smith }
900*d763cef2SBarry Smith 
901*d763cef2SBarry Smith #undef __FUNC__
902*d763cef2SBarry Smith #define __FUNC__ "TSStep"
903*d763cef2SBarry Smith /*@
904*d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
905*d763cef2SBarry Smith 
906*d763cef2SBarry Smith    Collective on TS
907*d763cef2SBarry Smith 
908*d763cef2SBarry Smith    Input Parameter:
909*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
910*d763cef2SBarry Smith 
911*d763cef2SBarry Smith    Output Parameters:
912*d763cef2SBarry Smith +  steps - number of iterations until termination
913*d763cef2SBarry Smith -  time - time until termination
914*d763cef2SBarry Smith 
915*d763cef2SBarry Smith    Level: beginner
916*d763cef2SBarry Smith 
917*d763cef2SBarry Smith .keywords: TS, timestep, solve
918*d763cef2SBarry Smith 
919*d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
920*d763cef2SBarry Smith @*/
921*d763cef2SBarry Smith int TSStep(TS ts,int *steps,double *time)
922*d763cef2SBarry Smith {
923*d763cef2SBarry Smith   int ierr,flg;
924*d763cef2SBarry Smith 
925*d763cef2SBarry Smith   PetscFunctionBegin;
926*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
927*d763cef2SBarry Smith   if (!ts->setupcalled) {ierr = TSSetUp(ts); CHKERRQ(ierr);}
928*d763cef2SBarry Smith   PLogEventBegin(TS_Step,ts,0,0,0);
929*d763cef2SBarry Smith   ierr = (*(ts)->step)(ts,steps,time); CHKERRQ(ierr);
930*d763cef2SBarry Smith   PLogEventEnd(TS_Step,ts,0,0,0);
931*d763cef2SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-ts_view",&flg); CHKERRQ(ierr);
932*d763cef2SBarry Smith   if (flg) {
933*d763cef2SBarry Smith     ierr = TSView(ts,VIEWER_STDOUT_WORLD); CHKERRQ(ierr);
934*d763cef2SBarry Smith   }
935*d763cef2SBarry Smith   PetscFunctionReturn(0);
936*d763cef2SBarry Smith }
937*d763cef2SBarry Smith 
938*d763cef2SBarry Smith #undef __FUNC__
939*d763cef2SBarry Smith #define __FUNC__ "TSMonitor"
940*d763cef2SBarry Smith /*
941*d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
942*d763cef2SBarry Smith */
943*d763cef2SBarry Smith int TSMonitor(TS ts,int step,double time,Vec x)
944*d763cef2SBarry Smith {
945*d763cef2SBarry Smith   int i,ierr,n = ts->numbermonitors;
946*d763cef2SBarry Smith 
947*d763cef2SBarry Smith   PetscFunctionBegin;
948*d763cef2SBarry Smith   for ( i=0; i<n; i++ ) {
949*d763cef2SBarry Smith     ierr = (*ts->monitor[i])(ts,step,time,x,ts->monitorcontext[i]);CHKERRQ(ierr);
950*d763cef2SBarry Smith   }
951*d763cef2SBarry Smith   PetscFunctionReturn(0);
952*d763cef2SBarry Smith }
953*d763cef2SBarry Smith 
954*d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
955*d763cef2SBarry Smith 
956*d763cef2SBarry Smith /*@C
957*d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
958*d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
959*d763cef2SBarry Smith 
960*d763cef2SBarry Smith    Collective on TS
961*d763cef2SBarry Smith 
962*d763cef2SBarry Smith    Input Parameters:
963*d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
964*d763cef2SBarry Smith .  label - the title to put in the title bar
965*d763cef2SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of
966*d763cef2SBarry Smith           the window
967*d763cef2SBarry Smith -  m, n - the screen width and height in pixels
968*d763cef2SBarry Smith 
969*d763cef2SBarry Smith    Output Parameter:
970*d763cef2SBarry Smith .  draw - the drawing context
971*d763cef2SBarry Smith 
972*d763cef2SBarry Smith    Options Database Key:
973*d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
974*d763cef2SBarry Smith 
975*d763cef2SBarry Smith    Notes:
976*d763cef2SBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not DrawLGDestroy().
977*d763cef2SBarry Smith 
978*d763cef2SBarry Smith    Level: intermediate
979*d763cef2SBarry Smith 
980*d763cef2SBarry Smith .keywords: TS, monitor, line graph, residual, create
981*d763cef2SBarry Smith 
982*d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
983*d763cef2SBarry Smith @*/
984*d763cef2SBarry Smith int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,
985*d763cef2SBarry Smith                        int n, DrawLG *draw)
986*d763cef2SBarry Smith {
987*d763cef2SBarry Smith   Draw win;
988*d763cef2SBarry Smith   int  ierr;
989*d763cef2SBarry Smith 
990*d763cef2SBarry Smith   PetscFunctionBegin;
991*d763cef2SBarry Smith   ierr = DrawOpenX(PETSC_COMM_SELF,host,label,x,y,m,n,&win); CHKERRQ(ierr);
992*d763cef2SBarry Smith   ierr = DrawLGCreate(win,1,draw); CHKERRQ(ierr);
993*d763cef2SBarry Smith   ierr = DrawLGIndicateDataPoints(*draw); CHKERRQ(ierr);
994*d763cef2SBarry Smith 
995*d763cef2SBarry Smith   PLogObjectParent(*draw,win);
996*d763cef2SBarry Smith   PetscFunctionReturn(0);
997*d763cef2SBarry Smith }
998*d763cef2SBarry Smith 
999*d763cef2SBarry Smith #undef __FUNC__
1000*d763cef2SBarry Smith #define __FUNC__ "TSLGMonitor"
1001*d763cef2SBarry Smith int TSLGMonitor(TS ts,int n,double time,Vec v,void *monctx)
1002*d763cef2SBarry Smith {
1003*d763cef2SBarry Smith   DrawLG lg = (DrawLG) monctx;
1004*d763cef2SBarry Smith   double x,y = time;
1005*d763cef2SBarry Smith   int    ierr;
1006*d763cef2SBarry Smith 
1007*d763cef2SBarry Smith   PetscFunctionBegin;
1008*d763cef2SBarry Smith   if (!n) {ierr = DrawLGReset(lg);CHKERRQ(ierr);}
1009*d763cef2SBarry Smith   x = (double) n;
1010*d763cef2SBarry Smith   ierr = DrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1011*d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1012*d763cef2SBarry Smith     ierr = DrawLGDraw(lg);CHKERRQ(ierr);
1013*d763cef2SBarry Smith   }
1014*d763cef2SBarry Smith   PetscFunctionReturn(0);
1015*d763cef2SBarry Smith }
1016*d763cef2SBarry Smith 
1017*d763cef2SBarry Smith #undef __FUNC__
1018*d763cef2SBarry Smith #define __FUNC__ "TSLGMonitorDestroy"
1019*d763cef2SBarry Smith /*@C
1020*d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1021*d763cef2SBarry Smith    with TSLGMonitorCreate().
1022*d763cef2SBarry Smith 
1023*d763cef2SBarry Smith    Collective on DrawLG
1024*d763cef2SBarry Smith 
1025*d763cef2SBarry Smith    Input Parameter:
1026*d763cef2SBarry Smith .  draw - the drawing context
1027*d763cef2SBarry Smith 
1028*d763cef2SBarry Smith    Level: intermediate
1029*d763cef2SBarry Smith 
1030*d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1031*d763cef2SBarry Smith 
1032*d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1033*d763cef2SBarry Smith @*/
1034*d763cef2SBarry Smith int TSLGMonitorDestroy(DrawLG drawlg)
1035*d763cef2SBarry Smith {
1036*d763cef2SBarry Smith   Draw draw;
1037*d763cef2SBarry Smith   int  ierr;
1038*d763cef2SBarry Smith 
1039*d763cef2SBarry Smith   PetscFunctionBegin;
1040*d763cef2SBarry Smith   ierr = DrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1041*d763cef2SBarry Smith   ierr = DrawDestroy(draw);CHKERRQ(ierr);
1042*d763cef2SBarry Smith   ierr = DrawLGDestroy(drawlg);CHKERRQ(ierr);
1043*d763cef2SBarry Smith   PetscFunctionReturn(0);
1044*d763cef2SBarry Smith }
1045*d763cef2SBarry Smith 
1046*d763cef2SBarry Smith #undef __FUNC__
1047*d763cef2SBarry Smith #define __FUNC__ "TSGetTime"
1048*d763cef2SBarry Smith /*@
1049*d763cef2SBarry Smith    TSGetTime - Gets the current time.
1050*d763cef2SBarry Smith 
1051*d763cef2SBarry Smith    Not Collective
1052*d763cef2SBarry Smith 
1053*d763cef2SBarry Smith    Input Parameter:
1054*d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1055*d763cef2SBarry Smith 
1056*d763cef2SBarry Smith    Output Parameter:
1057*d763cef2SBarry Smith .  t  - the current time
1058*d763cef2SBarry Smith 
1059*d763cef2SBarry Smith    Contributed by: Matthew Knepley
1060*d763cef2SBarry Smith 
1061*d763cef2SBarry Smith    Level: beginner
1062*d763cef2SBarry Smith 
1063*d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1064*d763cef2SBarry Smith 
1065*d763cef2SBarry Smith .keywords: TS, get, time
1066*d763cef2SBarry Smith @*/
1067*d763cef2SBarry Smith int TSGetTime(TS ts, double* t)
1068*d763cef2SBarry Smith {
1069*d763cef2SBarry Smith   PetscFunctionBegin;
1070*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1071*d763cef2SBarry Smith   *t = ts->ptime;
1072*d763cef2SBarry Smith   PetscFunctionReturn(0);
1073*d763cef2SBarry Smith }
1074*d763cef2SBarry Smith 
1075*d763cef2SBarry Smith #undef __FUNC__
1076*d763cef2SBarry Smith #define __FUNC__ "TSGetProblemType"
1077*d763cef2SBarry Smith /*@C
1078*d763cef2SBarry Smith    TSGetProblemType - Returns the problem type of a TS (timestepper) context.
1079*d763cef2SBarry Smith 
1080*d763cef2SBarry Smith    Not Collective
1081*d763cef2SBarry Smith 
1082*d763cef2SBarry Smith    Input Parameter:
1083*d763cef2SBarry Smith .  ts   - The TS context obtained from TSCreate()
1084*d763cef2SBarry Smith 
1085*d763cef2SBarry Smith    Output Parameter:
1086*d763cef2SBarry Smith .  type - The problem type, TS_LINEAR or TS_NONLINEAR
1087*d763cef2SBarry Smith 
1088*d763cef2SBarry Smith    Level: intermediate
1089*d763cef2SBarry Smith 
1090*d763cef2SBarry Smith    Contributed by: Matthew Knepley
1091*d763cef2SBarry Smith 
1092*d763cef2SBarry Smith .keywords: ts, get, type
1093*d763cef2SBarry Smith 
1094*d763cef2SBarry Smith @*/
1095*d763cef2SBarry Smith int TSGetProblemType(TS ts, TSProblemType *type)
1096*d763cef2SBarry Smith {
1097*d763cef2SBarry Smith   PetscFunctionBegin;
1098*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1099*d763cef2SBarry Smith   *type = ts->problem_type;
1100*d763cef2SBarry Smith   PetscFunctionReturn(0);
1101*d763cef2SBarry Smith }
1102*d763cef2SBarry Smith 
1103*d763cef2SBarry Smith #undef __FUNC__
1104*d763cef2SBarry Smith #define __FUNC__ "TSSetOptionsPrefix"
1105*d763cef2SBarry Smith /*@C
1106*d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1107*d763cef2SBarry Smith    TS options in the database.
1108*d763cef2SBarry Smith 
1109*d763cef2SBarry Smith    Collective on TS
1110*d763cef2SBarry Smith 
1111*d763cef2SBarry Smith    Input Parameter:
1112*d763cef2SBarry Smith +  ts     - The TS context
1113*d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1114*d763cef2SBarry Smith 
1115*d763cef2SBarry Smith    Notes:
1116*d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1117*d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1118*d763cef2SBarry Smith    hyphen.
1119*d763cef2SBarry Smith 
1120*d763cef2SBarry Smith    Contributed by: Matthew Knepley
1121*d763cef2SBarry Smith 
1122*d763cef2SBarry Smith    Level: advanced
1123*d763cef2SBarry Smith 
1124*d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1125*d763cef2SBarry Smith 
1126*d763cef2SBarry Smith .seealso: TSSetFromOptions()
1127*d763cef2SBarry Smith 
1128*d763cef2SBarry Smith @*/
1129*d763cef2SBarry Smith int TSSetOptionsPrefix(TS ts, char *prefix)
1130*d763cef2SBarry Smith {
1131*d763cef2SBarry Smith   int ierr;
1132*d763cef2SBarry Smith 
1133*d763cef2SBarry Smith   PetscFunctionBegin;
1134*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1135*d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject) ts, prefix); CHKERRQ(ierr);
1136*d763cef2SBarry Smith   switch(ts->problem_type) {
1137*d763cef2SBarry Smith     case TS_NONLINEAR:
1138*d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes, prefix);              CHKERRQ(ierr);
1139*d763cef2SBarry Smith       break;
1140*d763cef2SBarry Smith     case TS_LINEAR:
1141*d763cef2SBarry Smith       ierr = SLESSetOptionsPrefix(ts->sles, prefix);              CHKERRQ(ierr);
1142*d763cef2SBarry Smith       break;
1143*d763cef2SBarry Smith   }
1144*d763cef2SBarry Smith   PetscFunctionReturn(0);
1145*d763cef2SBarry Smith }
1146*d763cef2SBarry Smith 
1147*d763cef2SBarry Smith 
1148*d763cef2SBarry Smith #undef __FUNC__
1149*d763cef2SBarry Smith #define __FUNC__ "TSAppendOptionsPrefix"
1150*d763cef2SBarry Smith /*@C
1151*d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1152*d763cef2SBarry Smith    TS options in the database.
1153*d763cef2SBarry Smith 
1154*d763cef2SBarry Smith    Collective on TS
1155*d763cef2SBarry Smith 
1156*d763cef2SBarry Smith    Input Parameter:
1157*d763cef2SBarry Smith +  ts     - The TS context
1158*d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1159*d763cef2SBarry Smith 
1160*d763cef2SBarry Smith    Notes:
1161*d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1162*d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1163*d763cef2SBarry Smith    hyphen.
1164*d763cef2SBarry Smith 
1165*d763cef2SBarry Smith    Contributed by: Matthew Knepley
1166*d763cef2SBarry Smith 
1167*d763cef2SBarry Smith    Level: advanced
1168*d763cef2SBarry Smith 
1169*d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1170*d763cef2SBarry Smith 
1171*d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1172*d763cef2SBarry Smith 
1173*d763cef2SBarry Smith @*/
1174*d763cef2SBarry Smith int TSAppendOptionsPrefix(TS ts, char *prefix)
1175*d763cef2SBarry Smith {
1176*d763cef2SBarry Smith   int ierr;
1177*d763cef2SBarry Smith 
1178*d763cef2SBarry Smith   PetscFunctionBegin;
1179*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1180*d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject) ts, prefix); CHKERRQ(ierr);
1181*d763cef2SBarry Smith   switch(ts->problem_type) {
1182*d763cef2SBarry Smith     case TS_NONLINEAR:
1183*d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes, prefix);              CHKERRQ(ierr);
1184*d763cef2SBarry Smith       break;
1185*d763cef2SBarry Smith     case TS_LINEAR:
1186*d763cef2SBarry Smith       ierr = SLESAppendOptionsPrefix(ts->sles, prefix);              CHKERRQ(ierr);
1187*d763cef2SBarry Smith       break;
1188*d763cef2SBarry Smith   }
1189*d763cef2SBarry Smith   PetscFunctionReturn(0);
1190*d763cef2SBarry Smith }
1191*d763cef2SBarry Smith 
1192*d763cef2SBarry Smith #undef __FUNC__
1193*d763cef2SBarry Smith #define __FUNC__ "TSGetOptionsPrefix"
1194*d763cef2SBarry Smith /*@C
1195*d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1196*d763cef2SBarry Smith    TS options in the database.
1197*d763cef2SBarry Smith 
1198*d763cef2SBarry Smith    Not Collective
1199*d763cef2SBarry Smith 
1200*d763cef2SBarry Smith    Input Parameter:
1201*d763cef2SBarry Smith .  ts - The TS context
1202*d763cef2SBarry Smith 
1203*d763cef2SBarry Smith    Output Parameter:
1204*d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1205*d763cef2SBarry Smith 
1206*d763cef2SBarry Smith    Contributed by: Matthew Knepley
1207*d763cef2SBarry Smith 
1208*d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1209*d763cef2SBarry Smith    sufficient length to hold the prefix.
1210*d763cef2SBarry Smith 
1211*d763cef2SBarry Smith    Level: intermediate
1212*d763cef2SBarry Smith 
1213*d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1214*d763cef2SBarry Smith 
1215*d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1216*d763cef2SBarry Smith @*/
1217*d763cef2SBarry Smith int TSGetOptionsPrefix(TS ts, char **prefix)
1218*d763cef2SBarry Smith {
1219*d763cef2SBarry Smith   int ierr;
1220*d763cef2SBarry Smith 
1221*d763cef2SBarry Smith   PetscFunctionBegin;
1222*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1223*d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject) ts, prefix); CHKERRQ(ierr);
1224*d763cef2SBarry Smith   PetscFunctionReturn(0);
1225*d763cef2SBarry Smith }
1226*d763cef2SBarry Smith 
1227*d763cef2SBarry Smith #undef __FUNC__
1228*d763cef2SBarry Smith #define __FUNC__ "TSGetRHSMatrix"
1229*d763cef2SBarry Smith /*@C
1230*d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1231*d763cef2SBarry Smith 
1232*d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1233*d763cef2SBarry Smith 
1234*d763cef2SBarry Smith    Input Parameter:
1235*d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1236*d763cef2SBarry Smith 
1237*d763cef2SBarry Smith    Output Parameters:
1238*d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1239*d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1240*d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1241*d763cef2SBarry Smith 
1242*d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1243*d763cef2SBarry Smith 
1244*d763cef2SBarry Smith    Contributed by: Matthew Knepley
1245*d763cef2SBarry Smith 
1246*d763cef2SBarry Smith    Level: intermediate
1247*d763cef2SBarry Smith 
1248*d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1249*d763cef2SBarry Smith 
1250*d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1251*d763cef2SBarry Smith 
1252*d763cef2SBarry Smith @*/
1253*d763cef2SBarry Smith int TSGetRHSMatrix(TS ts, Mat *A, Mat *M, void **ctx)
1254*d763cef2SBarry Smith {
1255*d763cef2SBarry Smith   PetscFunctionBegin;
1256*d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1257*d763cef2SBarry Smith   if (A)   *A = ts->A;
1258*d763cef2SBarry Smith   if (M)   *M = ts->B;
1259*d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1260*d763cef2SBarry Smith   PetscFunctionReturn(0);
1261*d763cef2SBarry Smith }
1262*d763cef2SBarry Smith 
1263*d763cef2SBarry Smith #undef __FUNC__
1264*d763cef2SBarry Smith #define __FUNC__ "TSGetRHSJacobian"
1265*d763cef2SBarry Smith /*@C
1266*d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1267*d763cef2SBarry Smith 
1268*d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1269*d763cef2SBarry Smith 
1270*d763cef2SBarry Smith    Input Parameter:
1271*d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1272*d763cef2SBarry Smith 
1273*d763cef2SBarry Smith    Output Parameters:
1274*d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1275*d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1276*d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1277*d763cef2SBarry Smith 
1278*d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1279*d763cef2SBarry Smith 
1280*d763cef2SBarry Smith    Contributed by: Matthew Knepley
1281*d763cef2SBarry Smith 
1282*d763cef2SBarry Smith    Level: intermediate
1283*d763cef2SBarry Smith 
1284*d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1285*d763cef2SBarry Smith 
1286*d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1287*d763cef2SBarry Smith @*/
1288*d763cef2SBarry Smith int TSGetRHSJacobian(TS ts, Mat *J, Mat *M, void **ctx)
1289*d763cef2SBarry Smith {
1290*d763cef2SBarry Smith   int ierr;
1291*d763cef2SBarry Smith 
1292*d763cef2SBarry Smith   PetscFunctionBegin;
1293*d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts, J, M, ctx);CHKERRQ(ierr);
1294*d763cef2SBarry Smith   PetscFunctionReturn(0);
1295*d763cef2SBarry Smith }
1296*d763cef2SBarry Smith 
1297*d763cef2SBarry Smith /*MC
1298*d763cef2SBarry Smith    TSRegister - Adds a method to the timestepping solver package.
1299*d763cef2SBarry Smith 
1300*d763cef2SBarry Smith    Synopsis:
1301*d763cef2SBarry Smith 
1302*d763cef2SBarry Smith    TSRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(TS))
1303*d763cef2SBarry Smith 
1304*d763cef2SBarry Smith    Not collective
1305*d763cef2SBarry Smith 
1306*d763cef2SBarry Smith    Input Parameters:
1307*d763cef2SBarry Smith +  name_solver - name of a new user-defined solver
1308*d763cef2SBarry Smith .  path - path (either absolute or relative) the library containing this solver
1309*d763cef2SBarry Smith .  name_create - name of routine to create method context
1310*d763cef2SBarry Smith -  routine_create - routine to create method context
1311*d763cef2SBarry Smith 
1312*d763cef2SBarry Smith    Notes:
1313*d763cef2SBarry Smith    TSRegister() may be called multiple times to add several user-defined solvers.
1314*d763cef2SBarry Smith 
1315*d763cef2SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1316*d763cef2SBarry Smith    is ignored.
1317*d763cef2SBarry Smith 
1318*d763cef2SBarry Smith    Sample usage:
1319*d763cef2SBarry Smith .vb
1320*d763cef2SBarry Smith    TSRegister("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
1321*d763cef2SBarry Smith               "MySolverCreate",MySolverCreate);
1322*d763cef2SBarry Smith .ve
1323*d763cef2SBarry Smith 
1324*d763cef2SBarry Smith    Then, your solver can be chosen with the procedural interface via
1325*d763cef2SBarry Smith $     TSSetType(ts,"my_solver")
1326*d763cef2SBarry Smith    or at runtime via the option
1327*d763cef2SBarry Smith $     -ts_type my_solver
1328*d763cef2SBarry Smith 
1329*d763cef2SBarry Smith    Level: advanced
1330*d763cef2SBarry Smith 
1331*d763cef2SBarry Smith    $PETSC_ARCH and $BOPT occuring in pathname will be replaced with appropriate values.
1332*d763cef2SBarry Smith 
1333*d763cef2SBarry Smith .keywords: TS, register
1334*d763cef2SBarry Smith 
1335*d763cef2SBarry Smith .seealso: TSRegisterAll(), TSRegisterDestroy()
1336*d763cef2SBarry Smith M*/
1337*d763cef2SBarry Smith 
1338*d763cef2SBarry Smith #undef __FUNC__
1339*d763cef2SBarry Smith #define __FUNC__ "TSRegister_Private"
1340*d763cef2SBarry Smith int TSRegister_Private(char *sname,char *path,char *name,int (*function)(TS))
1341*d763cef2SBarry Smith {
1342*d763cef2SBarry Smith   char fullname[256];
1343*d763cef2SBarry Smith 
1344*d763cef2SBarry Smith   PetscFunctionBegin;
1345*d763cef2SBarry Smith   PetscStrcpy(fullname,path); PetscStrcat(fullname,":"); PetscStrcat(fullname,name);
1346*d763cef2SBarry Smith   FListAdd_Private(&TSList,sname,fullname,        (int (*)(void*))function);
1347*d763cef2SBarry Smith   PetscFunctionReturn(0);
1348*d763cef2SBarry Smith }
1349