xref: /petsc/src/ts/interface/ts.c (revision bdad233f16b1da2fc2320f4673f853efec1db8b6)
173f4d377SMatthew Knepley /* $Id: ts.c,v 1.43 2001/09/07 20:12:01 bsmith Exp $ */
2e090d566SSatish Balay #include "src/ts/tsimpl.h"        /*I "petscts.h"  I*/
3d763cef2SBarry Smith 
4d5ba7fb7SMatthew Knepley /* Logging support */
58ba1e511SMatthew Knepley int TS_COOKIE;
6d5ba7fb7SMatthew Knepley int TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
7d405a339SMatthew Knepley 
84a2ae208SSatish Balay #undef __FUNCT__
9*bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
10*bdad233fSMatthew Knepley /*
11*bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
12*bdad233fSMatthew Knepley 
13*bdad233fSMatthew Knepley   Collective on TS
14*bdad233fSMatthew Knepley 
15*bdad233fSMatthew Knepley   Input Parameter:
16*bdad233fSMatthew Knepley . ts - The ts
17*bdad233fSMatthew Knepley 
18*bdad233fSMatthew Knepley   Level: intermediate
19*bdad233fSMatthew Knepley 
20*bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
21*bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
22*bdad233fSMatthew Knepley */
23*bdad233fSMatthew Knepley static int TSSetTypeFromOptions(TS ts)
24*bdad233fSMatthew Knepley {
25*bdad233fSMatthew Knepley   PetscTruth opt;
26*bdad233fSMatthew Knepley   char      *defaultType;
27*bdad233fSMatthew Knepley   char       typeName[256];
28*bdad233fSMatthew Knepley   int        ierr;
29*bdad233fSMatthew Knepley 
30*bdad233fSMatthew Knepley   PetscFunctionBegin;
31*bdad233fSMatthew Knepley   if (ts->type_name != PETSC_NULL) {
32*bdad233fSMatthew Knepley     defaultType = ts->type_name;
33*bdad233fSMatthew Knepley   } else {
34*bdad233fSMatthew Knepley     defaultType = TS_EULER;
35*bdad233fSMatthew Knepley   }
36*bdad233fSMatthew Knepley 
37*bdad233fSMatthew Knepley   if (!TSRegisterAllCalled) {
38*bdad233fSMatthew Knepley     ierr = TSRegisterAll(PETSC_NULL);                                                                     CHKERRQ(ierr);
39*bdad233fSMatthew Knepley   }
40*bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
41*bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
42*bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);                                                                       CHKERRQ(ierr);
43*bdad233fSMatthew Knepley   } else {
44*bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);                                                                    CHKERRQ(ierr);
45*bdad233fSMatthew Knepley   }
46*bdad233fSMatthew Knepley   PetscFunctionReturn(0);
47*bdad233fSMatthew Knepley   return(0);
48*bdad233fSMatthew Knepley }
49*bdad233fSMatthew Knepley 
50*bdad233fSMatthew Knepley #undef __FUNCT__
51*bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
52*bdad233fSMatthew Knepley /*@
53*bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
54*bdad233fSMatthew Knepley 
55*bdad233fSMatthew Knepley    Collective on TS
56*bdad233fSMatthew Knepley 
57*bdad233fSMatthew Knepley    Input Parameter:
58*bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
59*bdad233fSMatthew Knepley 
60*bdad233fSMatthew Knepley    Options Database Keys:
61*bdad233fSMatthew Knepley +  -ts_type <type> - TS_EULER, TS_BEULER, TS_PVODE, TS_PSEUDO, TS_CRANK_NICHOLSON
62*bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
63*bdad233fSMatthew Knepley .  -ts_max_time time - maximum time to compute to
64*bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
65*bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
66*bdad233fSMatthew Knepley -  -ts_xmonitor - plot information at each timestep
67*bdad233fSMatthew Knepley 
68*bdad233fSMatthew Knepley    Level: beginner
69*bdad233fSMatthew Knepley 
70*bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
71*bdad233fSMatthew Knepley 
72*bdad233fSMatthew Knepley .seealso: TSGetType
73*bdad233fSMatthew Knepley @*/
74*bdad233fSMatthew Knepley int TSSetFromOptions(TS ts)
75*bdad233fSMatthew Knepley {
76*bdad233fSMatthew Knepley   PetscReal  dt;
77*bdad233fSMatthew Knepley   PetscTruth opt;
78*bdad233fSMatthew Knepley   int        ierr;
79*bdad233fSMatthew Knepley 
80*bdad233fSMatthew Knepley   PetscFunctionBegin;
81*bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
82*bdad233fSMatthew Knepley   ierr = PetscOptionsBegin(ts->comm, ts->prefix, "Time step options", "TS");                              CHKERRQ(ierr);
83*bdad233fSMatthew Knepley 
84*bdad233fSMatthew Knepley   /* Handle generic TS options */
85*bdad233fSMatthew Knepley   ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
86*bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
87*bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
88*bdad233fSMatthew Knepley   ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
89*bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
90*bdad233fSMatthew Knepley     ts->initial_time_step = ts->time_step = dt;
91*bdad233fSMatthew Knepley   }
92*bdad233fSMatthew Knepley 
93*bdad233fSMatthew Knepley   /* Monitor options */
94*bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_monitor","Monitor timestep size","TSDefaultMonitor",&opt);               CHKERRQ(ierr);
95*bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
96*bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSDefaultMonitor,PETSC_NULL,PETSC_NULL);                                     CHKERRQ(ierr);
97*bdad233fSMatthew Knepley     }
98*bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_xmonitor","Monitor timestep size graphically","TSLGMonitor",&opt);       CHKERRQ(ierr);
99*bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
100*bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSLGMonitor,PETSC_NULL,PETSC_NULL);                                          CHKERRQ(ierr);
101*bdad233fSMatthew Knepley     }
102*bdad233fSMatthew Knepley     ierr = PetscOptionsName("-ts_vecmonitor","Monitor solution graphically","TSVecViewMonitor",&opt);     CHKERRQ(ierr);
103*bdad233fSMatthew Knepley     if (opt == PETSC_TRUE) {
104*bdad233fSMatthew Knepley       ierr = TSSetMonitor(ts,TSVecViewMonitor,PETSC_NULL,PETSC_NULL);                                     CHKERRQ(ierr);
105*bdad233fSMatthew Knepley     }
106*bdad233fSMatthew Knepley 
107*bdad233fSMatthew Knepley   /* Handle TS type options */
108*bdad233fSMatthew Knepley   ierr = TSSetTypeFromOptions(ts);                                                                        CHKERRQ(ierr);
109*bdad233fSMatthew Knepley 
110*bdad233fSMatthew Knepley   /* Handle specific TS options */
111*bdad233fSMatthew Knepley   if (ts->ops->setfromoptions != PETSC_NULL) {
112*bdad233fSMatthew Knepley     ierr = (*ts->ops->setfromoptions)(ts);                                                                CHKERRQ(ierr);
113*bdad233fSMatthew Knepley   }
114*bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();                                                                               CHKERRQ(ierr);
115*bdad233fSMatthew Knepley 
116*bdad233fSMatthew Knepley   /* Handle subobject options */
117*bdad233fSMatthew Knepley   switch(ts->problem_type) {
118*bdad233fSMatthew Knepley   case TS_LINEAR:
119*bdad233fSMatthew Knepley     ierr = SLESSetFromOptions(ts->sles);                                                                  CHKERRQ(ierr);
120*bdad233fSMatthew Knepley     break;
121*bdad233fSMatthew Knepley   case TS_NONLINEAR:
122*bdad233fSMatthew Knepley     ierr = SNESSetFromOptions(ts->snes);                                                                  CHKERRQ(ierr);
123*bdad233fSMatthew Knepley     break;
124*bdad233fSMatthew Knepley   default:
125*bdad233fSMatthew Knepley     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", ts->problem_type);
126*bdad233fSMatthew Knepley   }
127*bdad233fSMatthew Knepley 
128*bdad233fSMatthew Knepley   ierr = TSViewFromOptions(ts, ts->name);                                                                 CHKERRQ(ierr);
129*bdad233fSMatthew Knepley   PetscFunctionReturn(0);
130*bdad233fSMatthew Knepley }
131*bdad233fSMatthew Knepley 
132*bdad233fSMatthew Knepley #undef  __FUNCT__
133*bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
134*bdad233fSMatthew Knepley /*@
135*bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
136*bdad233fSMatthew Knepley 
137*bdad233fSMatthew Knepley   Collective on TS
138*bdad233fSMatthew Knepley 
139*bdad233fSMatthew Knepley   Input Parameter:
140*bdad233fSMatthew Knepley . ts - The ts
141*bdad233fSMatthew Knepley 
142*bdad233fSMatthew Knepley   Level: intermediate
143*bdad233fSMatthew Knepley 
144*bdad233fSMatthew Knepley .keywords: TS, view, options, database
145*bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
146*bdad233fSMatthew Knepley @*/
147*bdad233fSMatthew Knepley int TSViewFromOptions(TS ts, char *title)
148*bdad233fSMatthew Knepley {
149*bdad233fSMatthew Knepley   PetscViewer viewer;
150*bdad233fSMatthew Knepley   PetscDraw   draw;
151*bdad233fSMatthew Knepley   PetscTruth  opt;
152*bdad233fSMatthew Knepley   char       *titleStr;
153*bdad233fSMatthew Knepley   char        typeName[1024];
154*bdad233fSMatthew Knepley   char        fileName[1024];
155*bdad233fSMatthew Knepley   int         len;
156*bdad233fSMatthew Knepley   int         ierr;
157*bdad233fSMatthew Knepley 
158*bdad233fSMatthew Knepley   PetscFunctionBegin;
159*bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);                                               CHKERRQ(ierr);
160*bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
161*bdad233fSMatthew Knepley     ierr = PetscOptionsGetString(ts->prefix, "-ts_view", typeName, 1024, &opt);                           CHKERRQ(ierr);
162*bdad233fSMatthew Knepley     ierr = PetscStrlen(typeName, &len);                                                                   CHKERRQ(ierr);
163*bdad233fSMatthew Knepley     if (len > 0) {
164*bdad233fSMatthew Knepley       ierr = PetscViewerCreate(ts->comm, &viewer);                                                        CHKERRQ(ierr);
165*bdad233fSMatthew Knepley       ierr = PetscViewerSetType(viewer, typeName);                                                        CHKERRQ(ierr);
166*bdad233fSMatthew Knepley       ierr = PetscOptionsGetString(ts->prefix, "-ts_view_file", fileName, 1024, &opt);                    CHKERRQ(ierr);
167*bdad233fSMatthew Knepley       if (opt == PETSC_TRUE) {
168*bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, fileName);                                                  CHKERRQ(ierr);
169*bdad233fSMatthew Knepley       } else {
170*bdad233fSMatthew Knepley         ierr = PetscViewerSetFilename(viewer, ts->name);                                                  CHKERRQ(ierr);
171*bdad233fSMatthew Knepley       }
172*bdad233fSMatthew Knepley       ierr = TSView(ts, viewer);                                                                          CHKERRQ(ierr);
173*bdad233fSMatthew Knepley       ierr = PetscViewerFlush(viewer);                                                                    CHKERRQ(ierr);
174*bdad233fSMatthew Knepley       ierr = PetscViewerDestroy(viewer);                                                                  CHKERRQ(ierr);
175*bdad233fSMatthew Knepley     } else {
176*bdad233fSMatthew Knepley       ierr = TSView(ts, PETSC_NULL);                                                                      CHKERRQ(ierr);
177*bdad233fSMatthew Knepley     }
178*bdad233fSMatthew Knepley   }
179*bdad233fSMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view_draw", &opt);                                          CHKERRQ(ierr);
180*bdad233fSMatthew Knepley   if (opt == PETSC_TRUE) {
181*bdad233fSMatthew Knepley     ierr = PetscViewerDrawOpen(ts->comm, 0, 0, 0, 0, 300, 300, &viewer);                                  CHKERRQ(ierr);
182*bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);                                                      CHKERRQ(ierr);
183*bdad233fSMatthew Knepley     if (title != PETSC_NULL) {
184*bdad233fSMatthew Knepley       titleStr = title;
185*bdad233fSMatthew Knepley     } else {
186*bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject) ts);                                                           CHKERRQ(ierr) ;
187*bdad233fSMatthew Knepley       titleStr = ts->name;
188*bdad233fSMatthew Knepley     }
189*bdad233fSMatthew Knepley     ierr = PetscDrawSetTitle(draw, titleStr);                                                             CHKERRQ(ierr);
190*bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);                                                                            CHKERRQ(ierr);
191*bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);                                                                      CHKERRQ(ierr);
192*bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);                                                                          CHKERRQ(ierr);
193*bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);                                                                    CHKERRQ(ierr);
194*bdad233fSMatthew Knepley   }
195*bdad233fSMatthew Knepley   PetscFunctionReturn(0);
196*bdad233fSMatthew Knepley }
197*bdad233fSMatthew Knepley 
198*bdad233fSMatthew Knepley #undef __FUNCT__
1994a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
200a7a1495cSBarry Smith /*@
2018c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
202a7a1495cSBarry Smith       set with TSSetRHSJacobian().
203a7a1495cSBarry Smith 
204a7a1495cSBarry Smith    Collective on TS and Vec
205a7a1495cSBarry Smith 
206a7a1495cSBarry Smith    Input Parameters:
207a7a1495cSBarry Smith +  ts - the SNES context
208a7a1495cSBarry Smith .  t - current timestep
209a7a1495cSBarry Smith -  x - input vector
210a7a1495cSBarry Smith 
211a7a1495cSBarry Smith    Output Parameters:
212a7a1495cSBarry Smith +  A - Jacobian matrix
213a7a1495cSBarry Smith .  B - optional preconditioning matrix
214a7a1495cSBarry Smith -  flag - flag indicating matrix structure
215a7a1495cSBarry Smith 
216a7a1495cSBarry Smith    Notes:
217a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
218a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
219a7a1495cSBarry Smith 
220a7a1495cSBarry Smith    See SLESSetOperators() for important information about setting the
221a7a1495cSBarry Smith    flag parameter.
222a7a1495cSBarry Smith 
223a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
224a7a1495cSBarry Smith 
225a7a1495cSBarry Smith    Level: developer
226a7a1495cSBarry Smith 
227a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
228a7a1495cSBarry Smith 
229a7a1495cSBarry Smith .seealso:  TSSetRHSJacobian(), SLESSetOperators()
230a7a1495cSBarry Smith @*/
23187828ca2SBarry Smith int TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
232a7a1495cSBarry Smith {
233a7a1495cSBarry Smith   int ierr;
234a7a1495cSBarry Smith 
235a7a1495cSBarry Smith   PetscFunctionBegin;
236a7a1495cSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
237a7a1495cSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE);
238a7a1495cSBarry Smith   PetscCheckSameComm(ts,X);
239a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
24029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
241a7a1495cSBarry Smith   }
242000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
243d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
244a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
245a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
246000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
247a7a1495cSBarry Smith     PetscStackPop;
248d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
249a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
250a7a1495cSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE);
251a7a1495cSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE);
252ef66eb69SBarry Smith   } else {
253ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
254ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255ef66eb69SBarry Smith     if (*A != *B) {
256ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258ef66eb69SBarry Smith     }
259ef66eb69SBarry Smith   }
260a7a1495cSBarry Smith   PetscFunctionReturn(0);
261a7a1495cSBarry Smith }
262a7a1495cSBarry Smith 
2634a2ae208SSatish Balay #undef __FUNCT__
2644a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
265d763cef2SBarry Smith /*
266d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
267d763cef2SBarry Smith 
268d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
269d763cef2SBarry Smith    this routine applies the matrix.
270d763cef2SBarry Smith */
27187828ca2SBarry Smith int TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
272d763cef2SBarry Smith {
273d763cef2SBarry Smith   int ierr;
274d763cef2SBarry Smith 
275d763cef2SBarry Smith   PetscFunctionBegin;
276d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
2777c922b88SBarry Smith   PetscValidHeader(x);
2787c922b88SBarry Smith   PetscValidHeader(y);
279d763cef2SBarry Smith 
280d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
281000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
282d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
283000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
284d763cef2SBarry Smith     PetscStackPop;
2857c922b88SBarry Smith   } else {
286000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
287d763cef2SBarry Smith       MatStructure flg;
288d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
289000e7ae3SMatthew Knepley       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
290d763cef2SBarry Smith       PetscStackPop;
291d763cef2SBarry Smith     }
292d763cef2SBarry Smith     ierr = MatMult(ts->A,x,y);CHKERRQ(ierr);
2937c922b88SBarry Smith   }
294d763cef2SBarry Smith 
295d763cef2SBarry Smith   /* apply user-provided boundary conditions (only needed if these are time dependent) */
296d763cef2SBarry Smith   ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr);
297d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
298d763cef2SBarry Smith 
299d763cef2SBarry Smith   PetscFunctionReturn(0);
300d763cef2SBarry Smith }
301d763cef2SBarry Smith 
3024a2ae208SSatish Balay #undef __FUNCT__
3034a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
304d763cef2SBarry Smith /*@C
305d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
306d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
307d763cef2SBarry Smith 
308d763cef2SBarry Smith     Collective on TS
309d763cef2SBarry Smith 
310d763cef2SBarry Smith     Input Parameters:
311d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
312d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
313d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
314d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
315d763cef2SBarry Smith 
316d763cef2SBarry Smith     Calling sequence of func:
31787828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
318d763cef2SBarry Smith 
319d763cef2SBarry Smith +   t - current timestep
320d763cef2SBarry Smith .   u - input vector
321d763cef2SBarry Smith .   F - function vector
322d763cef2SBarry Smith -   ctx - [optional] user-defined function context
323d763cef2SBarry Smith 
324d763cef2SBarry Smith     Important:
325d763cef2SBarry Smith     The user MUST call either this routine or TSSetRHSMatrix().
326d763cef2SBarry Smith 
327d763cef2SBarry Smith     Level: beginner
328d763cef2SBarry Smith 
329d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
330d763cef2SBarry Smith 
331d763cef2SBarry Smith .seealso: TSSetRHSMatrix()
332d763cef2SBarry Smith @*/
33387828ca2SBarry Smith int TSSetRHSFunction(TS ts,int (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
334d763cef2SBarry Smith {
335d763cef2SBarry Smith   PetscFunctionBegin;
336d763cef2SBarry Smith 
337d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
338d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
33929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
340d763cef2SBarry Smith   }
341000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
342d763cef2SBarry Smith   ts->funP             = ctx;
343d763cef2SBarry Smith   PetscFunctionReturn(0);
344d763cef2SBarry Smith }
345d763cef2SBarry Smith 
3464a2ae208SSatish Balay #undef __FUNCT__
3474a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSMatrix"
348d763cef2SBarry Smith /*@C
349d763cef2SBarry Smith    TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U.
350d763cef2SBarry Smith    Also sets the location to store A.
351d763cef2SBarry Smith 
352d763cef2SBarry Smith    Collective on TS
353d763cef2SBarry Smith 
354d763cef2SBarry Smith    Input Parameters:
355d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
356d763cef2SBarry Smith .  A   - matrix
357d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
358d763cef2SBarry Smith .  f   - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
359d763cef2SBarry Smith          if A is not a function of t.
360d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
361d763cef2SBarry Smith           matrix evaluation routine (may be PETSC_NULL)
362d763cef2SBarry Smith 
363d763cef2SBarry Smith    Calling sequence of func:
36487828ca2SBarry Smith $     func (TS ts,PetscReal t,Mat *A,Mat *B,int *flag,void *ctx);
365d763cef2SBarry Smith 
366d763cef2SBarry Smith +  t - current timestep
367d763cef2SBarry Smith .  A - matrix A, where U_t = A(t) U
368d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
369d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
370d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
371d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
372d763cef2SBarry Smith 
373d763cef2SBarry Smith    Notes:
374d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
375d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
376d763cef2SBarry Smith 
377d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
378d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
379d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
380d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
381d763cef2SBarry Smith    throughout the global iterations.
382d763cef2SBarry Smith 
383d763cef2SBarry Smith    Important:
384d763cef2SBarry Smith    The user MUST call either this routine or TSSetRHSFunction().
385d763cef2SBarry Smith 
386d763cef2SBarry Smith    Level: beginner
387d763cef2SBarry Smith 
388d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, matrix
389d763cef2SBarry Smith 
390d763cef2SBarry Smith .seealso: TSSetRHSFunction()
391d763cef2SBarry Smith @*/
39287828ca2SBarry Smith int TSSetRHSMatrix(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),void *ctx)
393d763cef2SBarry Smith {
394d763cef2SBarry Smith   PetscFunctionBegin;
395d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
396184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
397184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
398184914b5SBarry Smith   PetscCheckSameComm(ts,A);
399184914b5SBarry Smith   PetscCheckSameComm(ts,B);
400d763cef2SBarry Smith   if (ts->problem_type == TS_NONLINEAR) {
40129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for nonlinear problems; use TSSetRHSJacobian()");
402d763cef2SBarry Smith   }
403d763cef2SBarry Smith 
404000e7ae3SMatthew Knepley   ts->ops->rhsmatrix = f;
405d763cef2SBarry Smith   ts->jacP           = ctx;
406d763cef2SBarry Smith   ts->A              = A;
407d763cef2SBarry Smith   ts->B              = B;
408d763cef2SBarry Smith 
409d763cef2SBarry Smith   PetscFunctionReturn(0);
410d763cef2SBarry Smith }
411d763cef2SBarry Smith 
4124a2ae208SSatish Balay #undef __FUNCT__
4134a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
414d763cef2SBarry Smith /*@C
415d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
416d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
417d763cef2SBarry Smith 
418d763cef2SBarry Smith    Collective on TS
419d763cef2SBarry Smith 
420d763cef2SBarry Smith    Input Parameters:
421d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
422d763cef2SBarry Smith .  A   - Jacobian matrix
423d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
424d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
425d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
426d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
427d763cef2SBarry Smith 
428d763cef2SBarry Smith    Calling sequence of func:
42987828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
430d763cef2SBarry Smith 
431d763cef2SBarry Smith +  t - current timestep
432d763cef2SBarry Smith .  u - input vector
433d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
434d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
435d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
436d763cef2SBarry Smith           structure (same as flag in SLESSetOperators())
437d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
438d763cef2SBarry Smith 
439d763cef2SBarry Smith    Notes:
440d763cef2SBarry Smith    See SLESSetOperators() for important information about setting the flag
441d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
442d763cef2SBarry Smith 
443d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
444d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
445d763cef2SBarry Smith    completely new new matrix structure (not just different matrix elements)
446d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
447d763cef2SBarry Smith    throughout the global iterations.
448d763cef2SBarry Smith 
449d763cef2SBarry Smith    Level: beginner
450d763cef2SBarry Smith 
451d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
452d763cef2SBarry Smith 
453d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
454d763cef2SBarry Smith           SNESDefaultComputeJacobianColor()
455d763cef2SBarry Smith 
456d763cef2SBarry Smith @*/
45787828ca2SBarry Smith int TSSetRHSJacobian(TS ts,Mat A,Mat B,int (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
458d763cef2SBarry Smith {
459d763cef2SBarry Smith   PetscFunctionBegin;
460d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
461184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
462184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
463184914b5SBarry Smith   PetscCheckSameComm(ts,A);
464184914b5SBarry Smith   PetscCheckSameComm(ts,B);
465d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
46629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetRHSMatrix()");
467d763cef2SBarry Smith   }
468d763cef2SBarry Smith 
469000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
470d763cef2SBarry Smith   ts->jacP             = ctx;
471d763cef2SBarry Smith   ts->A                = A;
472d763cef2SBarry Smith   ts->B                = B;
473d763cef2SBarry Smith   PetscFunctionReturn(0);
474d763cef2SBarry Smith }
475d763cef2SBarry Smith 
4764a2ae208SSatish Balay #undef __FUNCT__
4774a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSBoundaryConditions"
478d763cef2SBarry Smith /*
479d763cef2SBarry Smith    TSComputeRHSBoundaryConditions - Evaluates the boundary condition function.
480d763cef2SBarry Smith 
481d763cef2SBarry Smith    Note: If the user did not provide a function but merely a matrix,
482d763cef2SBarry Smith    this routine applies the matrix.
483d763cef2SBarry Smith */
48487828ca2SBarry Smith int TSComputeRHSBoundaryConditions(TS ts,PetscReal t,Vec x)
485d763cef2SBarry Smith {
486d763cef2SBarry Smith   int ierr;
487d763cef2SBarry Smith 
488d763cef2SBarry Smith   PetscFunctionBegin;
489d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
490d763cef2SBarry Smith   PetscValidHeader(x);
491184914b5SBarry Smith   PetscCheckSameComm(ts,x);
492d763cef2SBarry Smith 
493000e7ae3SMatthew Knepley   if (ts->ops->rhsbc) {
494d763cef2SBarry Smith     PetscStackPush("TS user boundary condition function");
495000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr);
496d763cef2SBarry Smith     PetscStackPop;
497d763cef2SBarry Smith     PetscFunctionReturn(0);
498d763cef2SBarry Smith   }
499d763cef2SBarry Smith 
500d763cef2SBarry Smith   PetscFunctionReturn(0);
501d763cef2SBarry Smith }
502d763cef2SBarry Smith 
5034a2ae208SSatish Balay #undef __FUNCT__
5044a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSBoundaryConditions"
505d763cef2SBarry Smith /*@C
506d763cef2SBarry Smith     TSSetRHSBoundaryConditions - Sets the routine for evaluating the function,
507d763cef2SBarry Smith     boundary conditions for the function F.
508d763cef2SBarry Smith 
509d763cef2SBarry Smith     Collective on TS
510d763cef2SBarry Smith 
511d763cef2SBarry Smith     Input Parameters:
512d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
513d763cef2SBarry Smith .   f - routine for evaluating the boundary condition function
514d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
515d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
516d763cef2SBarry Smith 
517d763cef2SBarry Smith     Calling sequence of func:
51887828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec F,void *ctx);
519d763cef2SBarry Smith 
520d763cef2SBarry Smith +   t - current timestep
521d763cef2SBarry Smith .   F - function vector
522d763cef2SBarry Smith -   ctx - [optional] user-defined function context
523d763cef2SBarry Smith 
524d763cef2SBarry Smith     Level: intermediate
525d763cef2SBarry Smith 
526d763cef2SBarry Smith .keywords: TS, timestep, set, boundary conditions, function
527d763cef2SBarry Smith @*/
52887828ca2SBarry Smith int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,PetscReal,Vec,void*),void *ctx)
529d763cef2SBarry Smith {
530d763cef2SBarry Smith   PetscFunctionBegin;
531d763cef2SBarry Smith 
532d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
533d763cef2SBarry Smith   if (ts->problem_type != TS_LINEAR) {
53429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For linear problems only");
535d763cef2SBarry Smith   }
536000e7ae3SMatthew Knepley   ts->ops->rhsbc = f;
537d763cef2SBarry Smith   ts->bcP        = ctx;
538d763cef2SBarry Smith   PetscFunctionReturn(0);
539d763cef2SBarry Smith }
540d763cef2SBarry Smith 
5414a2ae208SSatish Balay #undef __FUNCT__
5424a2ae208SSatish Balay #define __FUNCT__ "TSView"
5437e2c5f70SBarry Smith /*@C
544d763cef2SBarry Smith     TSView - Prints the TS data structure.
545d763cef2SBarry Smith 
5464c49b128SBarry Smith     Collective on TS
547d763cef2SBarry Smith 
548d763cef2SBarry Smith     Input Parameters:
549d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
550d763cef2SBarry Smith -   viewer - visualization context
551d763cef2SBarry Smith 
552d763cef2SBarry Smith     Options Database Key:
553d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
554d763cef2SBarry Smith 
555d763cef2SBarry Smith     Notes:
556d763cef2SBarry Smith     The available visualization contexts include
557b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
558b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
559d763cef2SBarry Smith          output where only the first processor opens
560d763cef2SBarry Smith          the file.  All other processors send their
561d763cef2SBarry Smith          data to the first processor to print.
562d763cef2SBarry Smith 
563d763cef2SBarry Smith     The user can open an alternative visualization context with
564b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
565d763cef2SBarry Smith 
566d763cef2SBarry Smith     Level: beginner
567d763cef2SBarry Smith 
568d763cef2SBarry Smith .keywords: TS, timestep, view
569d763cef2SBarry Smith 
570b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
571d763cef2SBarry Smith @*/
572b0a32e0cSBarry Smith int TSView(TS ts,PetscViewer viewer)
573d763cef2SBarry Smith {
574d763cef2SBarry Smith   int        ierr;
575454a90a3SBarry Smith   char       *type;
5766831982aSBarry Smith   PetscTruth isascii,isstring;
577d763cef2SBarry Smith 
578d763cef2SBarry Smith   PetscFunctionBegin;
579d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
580b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ts->comm);
581b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
5826831982aSBarry Smith   PetscCheckSameComm(ts,viewer);
583fd16b177SBarry Smith 
584b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
585b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
5860f5bd95cSBarry Smith   if (isascii) {
587b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
588454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
589454a90a3SBarry Smith     if (type) {
590b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
591184914b5SBarry Smith     } else {
592b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
593184914b5SBarry Smith     }
594000e7ae3SMatthew Knepley     if (ts->ops->view) {
595b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
596000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
597b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
598d763cef2SBarry Smith     }
599b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr);
600b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%g\n",ts->max_time);CHKERRQ(ierr);
601d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
602b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr);
603d763cef2SBarry Smith     }
604b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr);
6050f5bd95cSBarry Smith   } else if (isstring) {
606454a90a3SBarry Smith     ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr);
607b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
608d763cef2SBarry Smith   }
609b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
610d763cef2SBarry Smith   if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);}
611d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
612b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
613d763cef2SBarry Smith   PetscFunctionReturn(0);
614d763cef2SBarry Smith }
615d763cef2SBarry Smith 
616d763cef2SBarry Smith 
6174a2ae208SSatish Balay #undef __FUNCT__
6184a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
619d763cef2SBarry Smith /*@C
620d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
621d763cef2SBarry Smith    the timesteppers.
622d763cef2SBarry Smith 
623d763cef2SBarry Smith    Collective on TS
624d763cef2SBarry Smith 
625d763cef2SBarry Smith    Input Parameters:
626d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
627d763cef2SBarry Smith -  usrP - optional user context
628d763cef2SBarry Smith 
629d763cef2SBarry Smith    Level: intermediate
630d763cef2SBarry Smith 
631d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
632d763cef2SBarry Smith 
633d763cef2SBarry Smith .seealso: TSGetApplicationContext()
634d763cef2SBarry Smith @*/
635d763cef2SBarry Smith int TSSetApplicationContext(TS ts,void *usrP)
636d763cef2SBarry Smith {
637d763cef2SBarry Smith   PetscFunctionBegin;
638d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
639d763cef2SBarry Smith   ts->user = usrP;
640d763cef2SBarry Smith   PetscFunctionReturn(0);
641d763cef2SBarry Smith }
642d763cef2SBarry Smith 
6434a2ae208SSatish Balay #undef __FUNCT__
6444a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
645d763cef2SBarry Smith /*@C
646d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
647d763cef2SBarry Smith     timestepper.
648d763cef2SBarry Smith 
649d763cef2SBarry Smith     Not Collective
650d763cef2SBarry Smith 
651d763cef2SBarry Smith     Input Parameter:
652d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
653d763cef2SBarry Smith 
654d763cef2SBarry Smith     Output Parameter:
655d763cef2SBarry Smith .   usrP - user context
656d763cef2SBarry Smith 
657d763cef2SBarry Smith     Level: intermediate
658d763cef2SBarry Smith 
659d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
660d763cef2SBarry Smith 
661d763cef2SBarry Smith .seealso: TSSetApplicationContext()
662d763cef2SBarry Smith @*/
663d763cef2SBarry Smith int TSGetApplicationContext(TS ts,void **usrP)
664d763cef2SBarry Smith {
665d763cef2SBarry Smith   PetscFunctionBegin;
666d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
667d763cef2SBarry Smith   *usrP = ts->user;
668d763cef2SBarry Smith   PetscFunctionReturn(0);
669d763cef2SBarry Smith }
670d763cef2SBarry Smith 
6714a2ae208SSatish Balay #undef __FUNCT__
6724a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
673d763cef2SBarry Smith /*@
674d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
675d763cef2SBarry Smith 
676d763cef2SBarry Smith    Not Collective
677d763cef2SBarry Smith 
678d763cef2SBarry Smith    Input Parameter:
679d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
680d763cef2SBarry Smith 
681d763cef2SBarry Smith    Output Parameter:
682d763cef2SBarry Smith .  iter - number steps so far
683d763cef2SBarry Smith 
684d763cef2SBarry Smith    Level: intermediate
685d763cef2SBarry Smith 
686d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
687d763cef2SBarry Smith @*/
688d763cef2SBarry Smith int TSGetTimeStepNumber(TS ts,int* iter)
689d763cef2SBarry Smith {
690d763cef2SBarry Smith   PetscFunctionBegin;
691d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
692d763cef2SBarry Smith   *iter = ts->steps;
693d763cef2SBarry Smith   PetscFunctionReturn(0);
694d763cef2SBarry Smith }
695d763cef2SBarry Smith 
6964a2ae208SSatish Balay #undef __FUNCT__
6974a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
698d763cef2SBarry Smith /*@
699d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
700d763cef2SBarry Smith    as well as the initial time.
701d763cef2SBarry Smith 
702d763cef2SBarry Smith    Collective on TS
703d763cef2SBarry Smith 
704d763cef2SBarry Smith    Input Parameters:
705d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
706d763cef2SBarry Smith .  initial_time - the initial time
707d763cef2SBarry Smith -  time_step - the size of the timestep
708d763cef2SBarry Smith 
709d763cef2SBarry Smith    Level: intermediate
710d763cef2SBarry Smith 
711d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
712d763cef2SBarry Smith 
713d763cef2SBarry Smith .keywords: TS, set, initial, timestep
714d763cef2SBarry Smith @*/
71587828ca2SBarry Smith int TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
716d763cef2SBarry Smith {
717d763cef2SBarry Smith   PetscFunctionBegin;
718d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
719d763cef2SBarry Smith   ts->time_step         = time_step;
720d763cef2SBarry Smith   ts->initial_time_step = time_step;
721d763cef2SBarry Smith   ts->ptime             = initial_time;
722d763cef2SBarry Smith   PetscFunctionReturn(0);
723d763cef2SBarry Smith }
724d763cef2SBarry Smith 
7254a2ae208SSatish Balay #undef __FUNCT__
7264a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
727d763cef2SBarry Smith /*@
728d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
729d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
730d763cef2SBarry Smith 
731d763cef2SBarry Smith    Collective on TS
732d763cef2SBarry Smith 
733d763cef2SBarry Smith    Input Parameters:
734d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
735d763cef2SBarry Smith -  time_step - the size of the timestep
736d763cef2SBarry Smith 
737d763cef2SBarry Smith    Level: intermediate
738d763cef2SBarry Smith 
739d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
740d763cef2SBarry Smith 
741d763cef2SBarry Smith .keywords: TS, set, timestep
742d763cef2SBarry Smith @*/
74387828ca2SBarry Smith int TSSetTimeStep(TS ts,PetscReal time_step)
744d763cef2SBarry Smith {
745d763cef2SBarry Smith   PetscFunctionBegin;
746d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
747d763cef2SBarry Smith   ts->time_step = time_step;
748d763cef2SBarry Smith   PetscFunctionReturn(0);
749d763cef2SBarry Smith }
750d763cef2SBarry Smith 
7514a2ae208SSatish Balay #undef __FUNCT__
7524a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
753d763cef2SBarry Smith /*@
754d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
755d763cef2SBarry Smith 
756d763cef2SBarry Smith    Not Collective
757d763cef2SBarry Smith 
758d763cef2SBarry Smith    Input Parameter:
759d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
760d763cef2SBarry Smith 
761d763cef2SBarry Smith    Output Parameter:
762d763cef2SBarry Smith .  dt - the current timestep size
763d763cef2SBarry Smith 
764d763cef2SBarry Smith    Level: intermediate
765d763cef2SBarry Smith 
766d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
767d763cef2SBarry Smith 
768d763cef2SBarry Smith .keywords: TS, get, timestep
769d763cef2SBarry Smith @*/
77087828ca2SBarry Smith int TSGetTimeStep(TS ts,PetscReal* dt)
771d763cef2SBarry Smith {
772d763cef2SBarry Smith   PetscFunctionBegin;
773d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
774d763cef2SBarry Smith   *dt = ts->time_step;
775d763cef2SBarry Smith   PetscFunctionReturn(0);
776d763cef2SBarry Smith }
777d763cef2SBarry Smith 
7784a2ae208SSatish Balay #undef __FUNCT__
7794a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
780d763cef2SBarry Smith /*@C
781d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
782d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
783d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
784d763cef2SBarry Smith    the solution at the next timestep has been calculated.
785d763cef2SBarry Smith 
786d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
787d763cef2SBarry Smith 
788d763cef2SBarry Smith    Input Parameter:
789d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
790d763cef2SBarry Smith 
791d763cef2SBarry Smith    Output Parameter:
792d763cef2SBarry Smith .  v - the vector containing the solution
793d763cef2SBarry Smith 
794d763cef2SBarry Smith    Level: intermediate
795d763cef2SBarry Smith 
796d763cef2SBarry Smith .seealso: TSGetTimeStep()
797d763cef2SBarry Smith 
798d763cef2SBarry Smith .keywords: TS, timestep, get, solution
799d763cef2SBarry Smith @*/
800d763cef2SBarry Smith int TSGetSolution(TS ts,Vec *v)
801d763cef2SBarry Smith {
802d763cef2SBarry Smith   PetscFunctionBegin;
803d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
804d763cef2SBarry Smith   *v = ts->vec_sol_always;
805d763cef2SBarry Smith   PetscFunctionReturn(0);
806d763cef2SBarry Smith }
807d763cef2SBarry Smith 
808*bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
8094a2ae208SSatish Balay #undef __FUNCT__
810*bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
811d763cef2SBarry Smith /*@C
812*bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
813d763cef2SBarry Smith 
814*bdad233fSMatthew Knepley   Not collective
815d763cef2SBarry Smith 
816*bdad233fSMatthew Knepley   Input Parameters:
817*bdad233fSMatthew Knepley + ts   - The TS
818*bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
819d763cef2SBarry Smith .vb
820d763cef2SBarry Smith          U_t = A U
821d763cef2SBarry Smith          U_t = A(t) U
822d763cef2SBarry Smith          U_t = F(t,U)
823d763cef2SBarry Smith .ve
824d763cef2SBarry Smith 
825d763cef2SBarry Smith    Level: beginner
826d763cef2SBarry Smith 
827*bdad233fSMatthew Knepley .keywords: TS, problem type
828*bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
829d763cef2SBarry Smith @*/
830*bdad233fSMatthew Knepley int TSSetProblemType(TS ts, TSProblemType type) {
831d763cef2SBarry Smith   PetscFunctionBegin;
832*bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
833*bdad233fSMatthew Knepley   ts->problem_type = type;
834d763cef2SBarry Smith   PetscFunctionReturn(0);
835d763cef2SBarry Smith }
836d763cef2SBarry Smith 
837*bdad233fSMatthew Knepley #undef __FUNCT__
838*bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
839*bdad233fSMatthew Knepley /*@C
840*bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
841*bdad233fSMatthew Knepley 
842*bdad233fSMatthew Knepley   Not collective
843*bdad233fSMatthew Knepley 
844*bdad233fSMatthew Knepley   Input Parameter:
845*bdad233fSMatthew Knepley . ts   - The TS
846*bdad233fSMatthew Knepley 
847*bdad233fSMatthew Knepley   Output Parameter:
848*bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
849*bdad233fSMatthew Knepley .vb
850*bdad233fSMatthew Knepley          U_t = A U
851*bdad233fSMatthew Knepley          U_t = A(t) U
852*bdad233fSMatthew Knepley          U_t = F(t,U)
853*bdad233fSMatthew Knepley .ve
854*bdad233fSMatthew Knepley 
855*bdad233fSMatthew Knepley    Level: beginner
856*bdad233fSMatthew Knepley 
857*bdad233fSMatthew Knepley .keywords: TS, problem type
858*bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
859*bdad233fSMatthew Knepley @*/
860*bdad233fSMatthew Knepley int TSGetProblemType(TS ts, TSProblemType *type) {
861*bdad233fSMatthew Knepley   PetscFunctionBegin;
862*bdad233fSMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
863*bdad233fSMatthew Knepley   PetscValidPointer(type);
864*bdad233fSMatthew Knepley   *type = ts->problem_type;
865*bdad233fSMatthew Knepley   PetscFunctionReturn(0);
866*bdad233fSMatthew Knepley }
867d763cef2SBarry Smith 
8684a2ae208SSatish Balay #undef __FUNCT__
8694a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
870d763cef2SBarry Smith /*@
871d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
872d763cef2SBarry Smith    of a timestepper.
873d763cef2SBarry Smith 
874d763cef2SBarry Smith    Collective on TS
875d763cef2SBarry Smith 
876d763cef2SBarry Smith    Input Parameter:
877d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
878d763cef2SBarry Smith 
879d763cef2SBarry Smith    Notes:
880d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
881d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
882d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
883d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
884d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
885d763cef2SBarry Smith 
886d763cef2SBarry Smith    Level: advanced
887d763cef2SBarry Smith 
888d763cef2SBarry Smith .keywords: TS, timestep, setup
889d763cef2SBarry Smith 
890d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
891d763cef2SBarry Smith @*/
892d763cef2SBarry Smith int TSSetUp(TS ts)
893d763cef2SBarry Smith {
894d763cef2SBarry Smith   int ierr;
895d763cef2SBarry Smith 
896d763cef2SBarry Smith   PetscFunctionBegin;
897d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
89829bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
899d763cef2SBarry Smith   if (!ts->type_name) {
900d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
901d763cef2SBarry Smith   }
902000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
903d763cef2SBarry Smith   ts->setupcalled = 1;
904d763cef2SBarry Smith   PetscFunctionReturn(0);
905d763cef2SBarry Smith }
906d763cef2SBarry Smith 
9074a2ae208SSatish Balay #undef __FUNCT__
9084a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
909d763cef2SBarry Smith /*@C
910d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
911d763cef2SBarry Smith    with TSCreate().
912d763cef2SBarry Smith 
913d763cef2SBarry Smith    Collective on TS
914d763cef2SBarry Smith 
915d763cef2SBarry Smith    Input Parameter:
916d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
917d763cef2SBarry Smith 
918d763cef2SBarry Smith    Level: beginner
919d763cef2SBarry Smith 
920d763cef2SBarry Smith .keywords: TS, timestepper, destroy
921d763cef2SBarry Smith 
922d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
923d763cef2SBarry Smith @*/
924d763cef2SBarry Smith int TSDestroy(TS ts)
925d763cef2SBarry Smith {
926329f5518SBarry Smith   int ierr,i;
927d763cef2SBarry Smith 
928d763cef2SBarry Smith   PetscFunctionBegin;
929d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
930d763cef2SBarry Smith   if (--ts->refct > 0) PetscFunctionReturn(0);
931d763cef2SBarry Smith 
932be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
9330f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
934be0abb6dSBarry Smith 
935d763cef2SBarry Smith   if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);}
936d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
937000e7ae3SMatthew Knepley   ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);
938329f5518SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
939329f5518SBarry Smith     if (ts->mdestroy[i]) {
940329f5518SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
941329f5518SBarry Smith     }
942329f5518SBarry Smith   }
943b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)ts);
944d763cef2SBarry Smith   PetscHeaderDestroy((PetscObject)ts);
945d763cef2SBarry Smith   PetscFunctionReturn(0);
946d763cef2SBarry Smith }
947d763cef2SBarry Smith 
9484a2ae208SSatish Balay #undef __FUNCT__
9494a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
950d763cef2SBarry Smith /*@C
951d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
952d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
953d763cef2SBarry Smith 
954d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
955d763cef2SBarry Smith 
956d763cef2SBarry Smith    Input Parameter:
957d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
958d763cef2SBarry Smith 
959d763cef2SBarry Smith    Output Parameter:
960d763cef2SBarry Smith .  snes - the nonlinear solver context
961d763cef2SBarry Smith 
962d763cef2SBarry Smith    Notes:
963d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
964d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
965d763cef2SBarry Smith    SLES, KSP, and PC contexts as well.
966d763cef2SBarry Smith 
967d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
968d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
969d763cef2SBarry Smith 
970d763cef2SBarry Smith    Level: beginner
971d763cef2SBarry Smith 
972d763cef2SBarry Smith .keywords: timestep, get, SNES
973d763cef2SBarry Smith @*/
974d763cef2SBarry Smith int TSGetSNES(TS ts,SNES *snes)
975d763cef2SBarry Smith {
976d763cef2SBarry Smith   PetscFunctionBegin;
977d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
97829bbc08cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetSLES()");
979d763cef2SBarry Smith   *snes = ts->snes;
980d763cef2SBarry Smith   PetscFunctionReturn(0);
981d763cef2SBarry Smith }
982d763cef2SBarry Smith 
9834a2ae208SSatish Balay #undef __FUNCT__
9844a2ae208SSatish Balay #define __FUNCT__ "TSGetSLES"
985d763cef2SBarry Smith /*@C
986d763cef2SBarry Smith    TSGetSLES - Returns the SLES (linear solver) associated with
987d763cef2SBarry Smith    a TS (timestepper) context.
988d763cef2SBarry Smith 
989d763cef2SBarry Smith    Not Collective, but SLES is parallel if TS is parallel
990d763cef2SBarry Smith 
991d763cef2SBarry Smith    Input Parameter:
992d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
993d763cef2SBarry Smith 
994d763cef2SBarry Smith    Output Parameter:
995d763cef2SBarry Smith .  sles - the nonlinear solver context
996d763cef2SBarry Smith 
997d763cef2SBarry Smith    Notes:
998d763cef2SBarry Smith    The user can then directly manipulate the SLES context to set various
999d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1000d763cef2SBarry Smith    KSP and PC contexts as well.
1001d763cef2SBarry Smith 
1002d763cef2SBarry Smith    TSGetSLES() does not work for integrators that do not use SLES;
1003d763cef2SBarry Smith    in this case TSGetSLES() returns PETSC_NULL in sles.
1004d763cef2SBarry Smith 
1005d763cef2SBarry Smith    Level: beginner
1006d763cef2SBarry Smith 
1007d763cef2SBarry Smith .keywords: timestep, get, SLES
1008d763cef2SBarry Smith @*/
1009d763cef2SBarry Smith int TSGetSLES(TS ts,SLES *sles)
1010d763cef2SBarry Smith {
1011d763cef2SBarry Smith   PetscFunctionBegin;
1012d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
101329bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1014d763cef2SBarry Smith   *sles = ts->sles;
1015d763cef2SBarry Smith   PetscFunctionReturn(0);
1016d763cef2SBarry Smith }
1017d763cef2SBarry Smith 
1018d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1019d763cef2SBarry Smith 
10204a2ae208SSatish Balay #undef __FUNCT__
10214a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1022d763cef2SBarry Smith /*@
1023d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1024d763cef2SBarry Smith    maximum time for iteration.
1025d763cef2SBarry Smith 
1026d763cef2SBarry Smith    Collective on TS
1027d763cef2SBarry Smith 
1028d763cef2SBarry Smith    Input Parameters:
1029d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1030d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1031d763cef2SBarry Smith -  maxtime - final time to iterate to
1032d763cef2SBarry Smith 
1033d763cef2SBarry Smith    Options Database Keys:
1034d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1035d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1036d763cef2SBarry Smith 
1037d763cef2SBarry Smith    Notes:
1038d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1039d763cef2SBarry Smith 
1040d763cef2SBarry Smith    Level: intermediate
1041d763cef2SBarry Smith 
1042d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1043d763cef2SBarry Smith @*/
104487828ca2SBarry Smith int TSSetDuration(TS ts,int maxsteps,PetscReal maxtime)
1045d763cef2SBarry Smith {
1046d763cef2SBarry Smith   PetscFunctionBegin;
1047d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1048d763cef2SBarry Smith   ts->max_steps = maxsteps;
1049d763cef2SBarry Smith   ts->max_time  = maxtime;
1050d763cef2SBarry Smith   PetscFunctionReturn(0);
1051d763cef2SBarry Smith }
1052d763cef2SBarry Smith 
10534a2ae208SSatish Balay #undef __FUNCT__
10544a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1055d763cef2SBarry Smith /*@
1056d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1057d763cef2SBarry Smith    for use by the TS routines.
1058d763cef2SBarry Smith 
1059d763cef2SBarry Smith    Collective on TS and Vec
1060d763cef2SBarry Smith 
1061d763cef2SBarry Smith    Input Parameters:
1062d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1063d763cef2SBarry Smith -  x - the solution vector
1064d763cef2SBarry Smith 
1065d763cef2SBarry Smith    Level: beginner
1066d763cef2SBarry Smith 
1067d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1068d763cef2SBarry Smith @*/
1069d763cef2SBarry Smith int TSSetSolution(TS ts,Vec x)
1070d763cef2SBarry Smith {
1071d763cef2SBarry Smith   PetscFunctionBegin;
1072d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1073d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1074d763cef2SBarry Smith   PetscFunctionReturn(0);
1075d763cef2SBarry Smith }
1076d763cef2SBarry Smith 
1077e74ef692SMatthew Knepley #undef __FUNCT__
1078e74ef692SMatthew Knepley #define __FUNCT__ "TSSetRhsBC"
1079000e7ae3SMatthew Knepley /*@
1080000e7ae3SMatthew Knepley   TSSetRhsBC - Sets the function which applies boundary conditions
1081000e7ae3SMatthew Knepley   to the Rhs of each system.
1082000e7ae3SMatthew Knepley 
1083000e7ae3SMatthew Knepley   Collective on TS
1084000e7ae3SMatthew Knepley 
1085000e7ae3SMatthew Knepley   Input Parameters:
1086000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1087000e7ae3SMatthew Knepley - func - The function
1088000e7ae3SMatthew Knepley 
1089000e7ae3SMatthew Knepley   Calling sequence of func:
1090000e7ae3SMatthew Knepley . func (TS ts, Vec rhs, void *ctx);
1091000e7ae3SMatthew Knepley 
1092000e7ae3SMatthew Knepley + rhs - The current rhs vector
1093000e7ae3SMatthew Knepley - ctx - The user-context
1094000e7ae3SMatthew Knepley 
1095000e7ae3SMatthew Knepley   Level: intermediate
1096000e7ae3SMatthew Knepley 
1097000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1098000e7ae3SMatthew Knepley @*/
1099000e7ae3SMatthew Knepley int TSSetRhsBC(TS ts, int (*func)(TS, Vec, void *))
1100000e7ae3SMatthew Knepley {
1101000e7ae3SMatthew Knepley   PetscFunctionBegin;
1102000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1103000e7ae3SMatthew Knepley   ts->ops->applyrhsbc = func;
1104000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1105000e7ae3SMatthew Knepley }
1106000e7ae3SMatthew Knepley 
1107e74ef692SMatthew Knepley #undef __FUNCT__
1108e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultRhsBC"
1109000e7ae3SMatthew Knepley /*@
1110000e7ae3SMatthew Knepley   TSDefaultRhsBC - The default boundary condition function which does nothing.
1111000e7ae3SMatthew Knepley 
1112000e7ae3SMatthew Knepley   Collective on TS
1113000e7ae3SMatthew Knepley 
1114000e7ae3SMatthew Knepley   Input Parameters:
1115000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1116000e7ae3SMatthew Knepley . rhs - The Rhs
1117000e7ae3SMatthew Knepley - ctx - The user-context
1118000e7ae3SMatthew Knepley 
1119000e7ae3SMatthew Knepley   Level: developer
1120000e7ae3SMatthew Knepley 
1121000e7ae3SMatthew Knepley .keywords: TS, Rhs, boundary conditions
1122000e7ae3SMatthew Knepley @*/
1123000e7ae3SMatthew Knepley int TSDefaultRhsBC(TS ts,  Vec rhs, void *ctx)
1124000e7ae3SMatthew Knepley {
1125000e7ae3SMatthew Knepley   PetscFunctionBegin;
1126000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1127000e7ae3SMatthew Knepley }
1128000e7ae3SMatthew Knepley 
1129e74ef692SMatthew Knepley #undef __FUNCT__
1130e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSystemMatrixBC"
1131000e7ae3SMatthew Knepley /*@
1132000e7ae3SMatthew Knepley   TSSetSystemMatrixBC - Sets the function which applies boundary conditions
1133000e7ae3SMatthew Knepley   to the system matrix and preconditioner of each system.
1134000e7ae3SMatthew Knepley 
1135000e7ae3SMatthew Knepley   Collective on TS
1136000e7ae3SMatthew Knepley 
1137000e7ae3SMatthew Knepley   Input Parameters:
1138000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1139000e7ae3SMatthew Knepley - func - The function
1140000e7ae3SMatthew Knepley 
1141000e7ae3SMatthew Knepley   Calling sequence of func:
1142000e7ae3SMatthew Knepley . func (TS ts, Mat A, Mat B, void *ctx);
1143000e7ae3SMatthew Knepley 
1144000e7ae3SMatthew Knepley + A   - The current system matrix
1145000e7ae3SMatthew Knepley . B   - The current preconditioner
1146000e7ae3SMatthew Knepley - ctx - The user-context
1147000e7ae3SMatthew Knepley 
1148000e7ae3SMatthew Knepley   Level: intermediate
1149000e7ae3SMatthew Knepley 
1150000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1151000e7ae3SMatthew Knepley @*/
1152000e7ae3SMatthew Knepley int TSSetSystemMatrixBC(TS ts, int (*func)(TS, Mat, Mat, void *))
1153000e7ae3SMatthew Knepley {
1154000e7ae3SMatthew Knepley   PetscFunctionBegin;
1155000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1156000e7ae3SMatthew Knepley   ts->ops->applymatrixbc = func;
1157000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1158000e7ae3SMatthew Knepley }
1159000e7ae3SMatthew Knepley 
1160e74ef692SMatthew Knepley #undef __FUNCT__
1161e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSystemMatrixBC"
1162000e7ae3SMatthew Knepley /*@
1163000e7ae3SMatthew Knepley   TSDefaultSystemMatrixBC - The default boundary condition function which
1164000e7ae3SMatthew Knepley   does nothing.
1165000e7ae3SMatthew Knepley 
1166000e7ae3SMatthew Knepley   Collective on TS
1167000e7ae3SMatthew Knepley 
1168000e7ae3SMatthew Knepley   Input Parameters:
1169000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1170000e7ae3SMatthew Knepley . A   - The system matrix
1171000e7ae3SMatthew Knepley . B   - The preconditioner
1172000e7ae3SMatthew Knepley - ctx - The user-context
1173000e7ae3SMatthew Knepley 
1174000e7ae3SMatthew Knepley   Level: developer
1175000e7ae3SMatthew Knepley 
1176000e7ae3SMatthew Knepley .keywords: TS, System matrix, boundary conditions
1177000e7ae3SMatthew Knepley @*/
1178000e7ae3SMatthew Knepley int TSDefaultSystemMatrixBC(TS ts, Mat A, Mat B, void *ctx)
1179000e7ae3SMatthew Knepley {
1180000e7ae3SMatthew Knepley   PetscFunctionBegin;
1181000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1182000e7ae3SMatthew Knepley }
1183000e7ae3SMatthew Knepley 
1184e74ef692SMatthew Knepley #undef __FUNCT__
1185e74ef692SMatthew Knepley #define __FUNCT__ "TSSetSolutionBC"
1186000e7ae3SMatthew Knepley /*@
1187000e7ae3SMatthew Knepley   TSSetSolutionBC - Sets the function which applies boundary conditions
1188000e7ae3SMatthew Knepley   to the solution of each system. This is necessary in nonlinear systems
1189000e7ae3SMatthew Knepley   which time dependent boundary conditions.
1190000e7ae3SMatthew Knepley 
1191000e7ae3SMatthew Knepley   Collective on TS
1192000e7ae3SMatthew Knepley 
1193000e7ae3SMatthew Knepley   Input Parameters:
1194000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1195000e7ae3SMatthew Knepley - func - The function
1196000e7ae3SMatthew Knepley 
1197000e7ae3SMatthew Knepley   Calling sequence of func:
1198000e7ae3SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
1199000e7ae3SMatthew Knepley 
1200000e7ae3SMatthew Knepley + sol - The current solution vector
1201000e7ae3SMatthew Knepley - ctx - The user-context
1202000e7ae3SMatthew Knepley 
1203000e7ae3SMatthew Knepley   Level: intermediate
1204000e7ae3SMatthew Knepley 
1205000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1206000e7ae3SMatthew Knepley @*/
1207000e7ae3SMatthew Knepley int TSSetSolutionBC(TS ts, int (*func)(TS, Vec, void *))
1208000e7ae3SMatthew Knepley {
1209000e7ae3SMatthew Knepley   PetscFunctionBegin;
1210000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1211000e7ae3SMatthew Knepley   ts->ops->applysolbc = func;
1212000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1213000e7ae3SMatthew Knepley }
1214000e7ae3SMatthew Knepley 
1215e74ef692SMatthew Knepley #undef __FUNCT__
1216e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultSolutionBC"
1217000e7ae3SMatthew Knepley /*@
1218000e7ae3SMatthew Knepley   TSDefaultSolutionBC - The default boundary condition function which
1219000e7ae3SMatthew Knepley   does nothing.
1220000e7ae3SMatthew Knepley 
1221000e7ae3SMatthew Knepley   Collective on TS
1222000e7ae3SMatthew Knepley 
1223000e7ae3SMatthew Knepley   Input Parameters:
1224000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1225000e7ae3SMatthew Knepley . sol - The solution
1226000e7ae3SMatthew Knepley - ctx - The user-context
1227000e7ae3SMatthew Knepley 
1228000e7ae3SMatthew Knepley   Level: developer
1229000e7ae3SMatthew Knepley 
1230000e7ae3SMatthew Knepley .keywords: TS, solution, boundary conditions
1231000e7ae3SMatthew Knepley @*/
1232000e7ae3SMatthew Knepley int TSDefaultSolutionBC(TS ts, Vec sol, void *ctx)
1233000e7ae3SMatthew Knepley {
1234000e7ae3SMatthew Knepley   PetscFunctionBegin;
1235000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1236000e7ae3SMatthew Knepley }
1237000e7ae3SMatthew Knepley 
1238e74ef692SMatthew Knepley #undef __FUNCT__
1239e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1240000e7ae3SMatthew Knepley /*@
1241000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1242000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1243000e7ae3SMatthew Knepley 
1244000e7ae3SMatthew Knepley   Collective on TS
1245000e7ae3SMatthew Knepley 
1246000e7ae3SMatthew Knepley   Input Parameters:
1247000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1248000e7ae3SMatthew Knepley - func - The function
1249000e7ae3SMatthew Knepley 
1250000e7ae3SMatthew Knepley   Calling sequence of func:
1251000e7ae3SMatthew Knepley . func (TS ts);
1252000e7ae3SMatthew Knepley 
1253000e7ae3SMatthew Knepley   Level: intermediate
1254000e7ae3SMatthew Knepley 
1255000e7ae3SMatthew Knepley .keywords: TS, timestep
1256000e7ae3SMatthew Knepley @*/
1257000e7ae3SMatthew Knepley int TSSetPreStep(TS ts, int (*func)(TS))
1258000e7ae3SMatthew Knepley {
1259000e7ae3SMatthew Knepley   PetscFunctionBegin;
1260000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1261000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1262000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1263000e7ae3SMatthew Knepley }
1264000e7ae3SMatthew Knepley 
1265e74ef692SMatthew Knepley #undef __FUNCT__
1266e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1267000e7ae3SMatthew Knepley /*@
1268000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1269000e7ae3SMatthew Knepley 
1270000e7ae3SMatthew Knepley   Collective on TS
1271000e7ae3SMatthew Knepley 
1272000e7ae3SMatthew Knepley   Input Parameters:
1273000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1274000e7ae3SMatthew Knepley 
1275000e7ae3SMatthew Knepley   Level: developer
1276000e7ae3SMatthew Knepley 
1277000e7ae3SMatthew Knepley .keywords: TS, timestep
1278000e7ae3SMatthew Knepley @*/
1279000e7ae3SMatthew Knepley int TSDefaultPreStep(TS ts)
1280000e7ae3SMatthew Knepley {
1281000e7ae3SMatthew Knepley   PetscFunctionBegin;
1282000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1283000e7ae3SMatthew Knepley }
1284000e7ae3SMatthew Knepley 
1285e74ef692SMatthew Knepley #undef __FUNCT__
1286e74ef692SMatthew Knepley #define __FUNCT__ "TSSetUpdate"
1287000e7ae3SMatthew Knepley /*@
1288000e7ae3SMatthew Knepley   TSSetUpdate - Sets the general-purpose update function called
1289000e7ae3SMatthew Knepley   at the beginning of every time step. This function can change
1290000e7ae3SMatthew Knepley   the time step.
1291000e7ae3SMatthew Knepley 
1292000e7ae3SMatthew Knepley   Collective on TS
1293000e7ae3SMatthew Knepley 
1294000e7ae3SMatthew Knepley   Input Parameters:
1295000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1296000e7ae3SMatthew Knepley - func - The function
1297000e7ae3SMatthew Knepley 
1298000e7ae3SMatthew Knepley   Calling sequence of func:
1299000e7ae3SMatthew Knepley . func (TS ts, double t, double *dt);
1300000e7ae3SMatthew Knepley 
1301000e7ae3SMatthew Knepley + t   - The current time
1302000e7ae3SMatthew Knepley - dt  - The current time step
1303000e7ae3SMatthew Knepley 
1304000e7ae3SMatthew Knepley   Level: intermediate
1305000e7ae3SMatthew Knepley 
1306000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1307000e7ae3SMatthew Knepley @*/
1308000e7ae3SMatthew Knepley int TSSetUpdate(TS ts, int (*func)(TS, double, double *))
1309000e7ae3SMatthew Knepley {
1310000e7ae3SMatthew Knepley   PetscFunctionBegin;
1311000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1312000e7ae3SMatthew Knepley   ts->ops->update = func;
1313000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1314000e7ae3SMatthew Knepley }
1315000e7ae3SMatthew Knepley 
1316e74ef692SMatthew Knepley #undef __FUNCT__
1317e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultUpdate"
1318000e7ae3SMatthew Knepley /*@
1319000e7ae3SMatthew Knepley   TSDefaultUpdate - The default update function which does nothing.
1320000e7ae3SMatthew Knepley 
1321000e7ae3SMatthew Knepley   Collective on TS
1322000e7ae3SMatthew Knepley 
1323000e7ae3SMatthew Knepley   Input Parameters:
1324000e7ae3SMatthew Knepley + ts  - The TS context obtained from TSCreate()
1325000e7ae3SMatthew Knepley - t   - The current time
1326000e7ae3SMatthew Knepley 
1327000e7ae3SMatthew Knepley   Output Parameters:
1328000e7ae3SMatthew Knepley . dt  - The current time step
1329000e7ae3SMatthew Knepley 
1330000e7ae3SMatthew Knepley   Level: developer
1331000e7ae3SMatthew Knepley 
1332000e7ae3SMatthew Knepley .keywords: TS, update, timestep
1333000e7ae3SMatthew Knepley @*/
1334000e7ae3SMatthew Knepley int TSDefaultUpdate(TS ts, double t, double *dt)
1335000e7ae3SMatthew Knepley {
1336000e7ae3SMatthew Knepley   PetscFunctionBegin;
1337000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1338000e7ae3SMatthew Knepley }
1339000e7ae3SMatthew Knepley 
1340e74ef692SMatthew Knepley #undef __FUNCT__
1341e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1342000e7ae3SMatthew Knepley /*@
1343000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1344000e7ae3SMatthew Knepley   called once at the end of time stepping.
1345000e7ae3SMatthew Knepley 
1346000e7ae3SMatthew Knepley   Collective on TS
1347000e7ae3SMatthew Knepley 
1348000e7ae3SMatthew Knepley   Input Parameters:
1349000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1350000e7ae3SMatthew Knepley - func - The function
1351000e7ae3SMatthew Knepley 
1352000e7ae3SMatthew Knepley   Calling sequence of func:
1353000e7ae3SMatthew Knepley . func (TS ts);
1354000e7ae3SMatthew Knepley 
1355000e7ae3SMatthew Knepley   Level: intermediate
1356000e7ae3SMatthew Knepley 
1357000e7ae3SMatthew Knepley .keywords: TS, timestep
1358000e7ae3SMatthew Knepley @*/
1359000e7ae3SMatthew Knepley int TSSetPostStep(TS ts, int (*func)(TS))
1360000e7ae3SMatthew Knepley {
1361000e7ae3SMatthew Knepley   PetscFunctionBegin;
1362000e7ae3SMatthew Knepley   PetscValidHeaderSpecific(ts, TS_COOKIE);
1363000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1364000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1365000e7ae3SMatthew Knepley }
1366000e7ae3SMatthew Knepley 
1367e74ef692SMatthew Knepley #undef __FUNCT__
1368e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1369000e7ae3SMatthew Knepley /*@
1370000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1371000e7ae3SMatthew Knepley 
1372000e7ae3SMatthew Knepley   Collective on TS
1373000e7ae3SMatthew Knepley 
1374000e7ae3SMatthew Knepley   Input Parameters:
1375000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1376000e7ae3SMatthew Knepley 
1377000e7ae3SMatthew Knepley   Level: developer
1378000e7ae3SMatthew Knepley 
1379000e7ae3SMatthew Knepley .keywords: TS, timestep
1380000e7ae3SMatthew Knepley @*/
1381000e7ae3SMatthew Knepley int TSDefaultPostStep(TS ts)
1382000e7ae3SMatthew Knepley {
1383000e7ae3SMatthew Knepley   PetscFunctionBegin;
1384000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1385000e7ae3SMatthew Knepley }
1386000e7ae3SMatthew Knepley 
1387d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1388d763cef2SBarry Smith 
13894a2ae208SSatish Balay #undef __FUNCT__
13904a2ae208SSatish Balay #define __FUNCT__ "TSSetMonitor"
1391d763cef2SBarry Smith /*@C
1392d763cef2SBarry Smith    TSSetMonitor - Sets an ADDITIONAL function that is to be used at every
1393d763cef2SBarry Smith    timestep to display the iteration's  progress.
1394d763cef2SBarry Smith 
1395d763cef2SBarry Smith    Collective on TS
1396d763cef2SBarry Smith 
1397d763cef2SBarry Smith    Input Parameters:
1398d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1399d763cef2SBarry Smith .  func - monitoring routine
1400329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1401b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1402b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1403b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1404d763cef2SBarry Smith 
1405d763cef2SBarry Smith    Calling sequence of func:
140687828ca2SBarry Smith $    int func(TS ts,int steps,PetscReal time,Vec x,void *mctx)
1407d763cef2SBarry Smith 
1408d763cef2SBarry Smith +    ts - the TS context
1409d763cef2SBarry Smith .    steps - iteration number
14101f06c33eSBarry Smith .    time - current time
1411d763cef2SBarry Smith .    x - current iterate
1412d763cef2SBarry Smith -    mctx - [optional] monitoring context
1413d763cef2SBarry Smith 
1414d763cef2SBarry Smith    Notes:
1415d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1416d763cef2SBarry Smith    already has been loaded.
1417d763cef2SBarry Smith 
1418d763cef2SBarry Smith    Level: intermediate
1419d763cef2SBarry Smith 
1420d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1421d763cef2SBarry Smith 
1422d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSClearMonitor()
1423d763cef2SBarry Smith @*/
142487828ca2SBarry Smith int TSSetMonitor(TS ts,int (*monitor)(TS,int,PetscReal,Vec,void*),void *mctx,int (*mdestroy)(void*))
1425d763cef2SBarry Smith {
1426d763cef2SBarry Smith   PetscFunctionBegin;
1427d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1428d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
142929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1430d763cef2SBarry Smith   }
1431d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1432329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1433d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1434d763cef2SBarry Smith   PetscFunctionReturn(0);
1435d763cef2SBarry Smith }
1436d763cef2SBarry Smith 
14374a2ae208SSatish Balay #undef __FUNCT__
14384a2ae208SSatish Balay #define __FUNCT__ "TSClearMonitor"
1439d763cef2SBarry Smith /*@C
1440d763cef2SBarry Smith    TSClearMonitor - Clears all the monitors that have been set on a time-step object.
1441d763cef2SBarry Smith 
1442d763cef2SBarry Smith    Collective on TS
1443d763cef2SBarry Smith 
1444d763cef2SBarry Smith    Input Parameters:
1445d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1446d763cef2SBarry Smith 
1447d763cef2SBarry Smith    Notes:
1448d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1449d763cef2SBarry Smith 
1450d763cef2SBarry Smith    Level: intermediate
1451d763cef2SBarry Smith 
1452d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1453d763cef2SBarry Smith 
1454d763cef2SBarry Smith .seealso: TSDefaultMonitor(), TSSetMonitor()
1455d763cef2SBarry Smith @*/
1456d763cef2SBarry Smith int TSClearMonitor(TS ts)
1457d763cef2SBarry Smith {
1458d763cef2SBarry Smith   PetscFunctionBegin;
1459d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1460d763cef2SBarry Smith   ts->numbermonitors = 0;
1461d763cef2SBarry Smith   PetscFunctionReturn(0);
1462d763cef2SBarry Smith }
1463d763cef2SBarry Smith 
14644a2ae208SSatish Balay #undef __FUNCT__
14654a2ae208SSatish Balay #define __FUNCT__ "TSDefaultMonitor"
146687828ca2SBarry Smith int TSDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx)
1467d763cef2SBarry Smith {
1468d132466eSBarry Smith   int ierr;
1469d132466eSBarry Smith 
1470d763cef2SBarry Smith   PetscFunctionBegin;
1471142b95e3SSatish Balay   ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1472d763cef2SBarry Smith   PetscFunctionReturn(0);
1473d763cef2SBarry Smith }
1474d763cef2SBarry Smith 
14754a2ae208SSatish Balay #undef __FUNCT__
14764a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1477d763cef2SBarry Smith /*@
1478d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1479d763cef2SBarry Smith 
1480d763cef2SBarry Smith    Collective on TS
1481d763cef2SBarry Smith 
1482d763cef2SBarry Smith    Input Parameter:
1483d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1484d763cef2SBarry Smith 
1485d763cef2SBarry Smith    Output Parameters:
1486d763cef2SBarry Smith +  steps - number of iterations until termination
1487142b95e3SSatish Balay -  ptime - time until termination
1488d763cef2SBarry Smith 
1489d763cef2SBarry Smith    Level: beginner
1490d763cef2SBarry Smith 
1491d763cef2SBarry Smith .keywords: TS, timestep, solve
1492d763cef2SBarry Smith 
1493d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1494d763cef2SBarry Smith @*/
149587828ca2SBarry Smith int TSStep(TS ts,int *steps,PetscReal *ptime)
1496d763cef2SBarry Smith {
1497d405a339SMatthew Knepley   PetscTruth opt;
1498f1af5d2fSBarry Smith   int        ierr;
1499d763cef2SBarry Smith 
1500d763cef2SBarry Smith   PetscFunctionBegin;
1501d763cef2SBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE);
1502d405a339SMatthew Knepley   if (!ts->setupcalled) {
1503d405a339SMatthew Knepley     ierr = TSSetUp(ts);                                                                                   CHKERRQ(ierr);
1504d405a339SMatthew Knepley   }
1505d405a339SMatthew Knepley 
1506d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);                                                        CHKERRQ(ierr);
1507d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);                                                                         CHKERRQ(ierr);
1508000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);                                                              CHKERRQ(ierr);
1509d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);                                                                        CHKERRQ(ierr);
1510d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);                                                          CHKERRQ(ierr);
1511d405a339SMatthew Knepley 
1512d405a339SMatthew Knepley   ierr = PetscOptionsHasName(ts->prefix, "-ts_view", &opt);                                               CHKERRQ(ierr);
1513d405a339SMatthew Knepley   if ((opt == PETSC_TRUE) && !PetscPreLoadingOn) {
1514d405a339SMatthew Knepley     ierr = TSView(ts, PETSC_VIEWER_STDOUT_WORLD);                                                         CHKERRQ(ierr);
1515d405a339SMatthew Knepley   }
1516d763cef2SBarry Smith   PetscFunctionReturn(0);
1517d763cef2SBarry Smith }
1518d763cef2SBarry Smith 
15194a2ae208SSatish Balay #undef __FUNCT__
15204a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1521d763cef2SBarry Smith /*
1522d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1523d763cef2SBarry Smith */
152487828ca2SBarry Smith int TSMonitor(TS ts,int step,PetscReal ptime,Vec x)
1525d763cef2SBarry Smith {
1526d763cef2SBarry Smith   int i,ierr,n = ts->numbermonitors;
1527d763cef2SBarry Smith 
1528d763cef2SBarry Smith   PetscFunctionBegin;
1529d763cef2SBarry Smith   for (i=0; i<n; i++) {
1530142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1531d763cef2SBarry Smith   }
1532d763cef2SBarry Smith   PetscFunctionReturn(0);
1533d763cef2SBarry Smith }
1534d763cef2SBarry Smith 
1535d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1536d763cef2SBarry Smith 
15374a2ae208SSatish Balay #undef __FUNCT__
15384a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorCreate"
1539d763cef2SBarry Smith /*@C
1540d763cef2SBarry Smith    TSLGMonitorCreate - Creates a line graph context for use with
1541d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1542d763cef2SBarry Smith 
1543d763cef2SBarry Smith    Collective on TS
1544d763cef2SBarry Smith 
1545d763cef2SBarry Smith    Input Parameters:
1546d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1547d763cef2SBarry Smith .  label - the title to put in the title bar
15487c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1549d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1550d763cef2SBarry Smith 
1551d763cef2SBarry Smith    Output Parameter:
1552d763cef2SBarry Smith .  draw - the drawing context
1553d763cef2SBarry Smith 
1554d763cef2SBarry Smith    Options Database Key:
1555d763cef2SBarry Smith .  -ts_xmonitor - automatically sets line graph monitor
1556d763cef2SBarry Smith 
1557d763cef2SBarry Smith    Notes:
1558b0a32e0cSBarry Smith    Use TSLGMonitorDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1559d763cef2SBarry Smith 
1560d763cef2SBarry Smith    Level: intermediate
1561d763cef2SBarry Smith 
15627c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1563d763cef2SBarry Smith 
1564d763cef2SBarry Smith .seealso: TSLGMonitorDestroy(), TSSetMonitor()
15657c922b88SBarry Smith 
1566d763cef2SBarry Smith @*/
1567b0a32e0cSBarry Smith int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1568d763cef2SBarry Smith {
1569b0a32e0cSBarry Smith   PetscDraw win;
1570d763cef2SBarry Smith   int       ierr;
1571d763cef2SBarry Smith 
1572d763cef2SBarry Smith   PetscFunctionBegin;
1573b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1574b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1575b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1576b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1577d763cef2SBarry Smith 
1578b0a32e0cSBarry Smith   PetscLogObjectParent(*draw,win);
1579d763cef2SBarry Smith   PetscFunctionReturn(0);
1580d763cef2SBarry Smith }
1581d763cef2SBarry Smith 
15824a2ae208SSatish Balay #undef __FUNCT__
15834a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitor"
158487828ca2SBarry Smith int TSLGMonitor(TS ts,int n,PetscReal ptime,Vec v,void *monctx)
1585d763cef2SBarry Smith {
1586b0a32e0cSBarry Smith   PetscDrawLG lg = (PetscDrawLG) monctx;
158787828ca2SBarry Smith   PetscReal      x,y = ptime;
1588d763cef2SBarry Smith   int         ierr;
1589d763cef2SBarry Smith 
1590d763cef2SBarry Smith   PetscFunctionBegin;
15917c922b88SBarry Smith   if (!monctx) {
15927c922b88SBarry Smith     MPI_Comm    comm;
1593b0a32e0cSBarry Smith     PetscViewer viewer;
15947c922b88SBarry Smith 
15957c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1596b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1597b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
15987c922b88SBarry Smith   }
15997c922b88SBarry Smith 
1600b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
160187828ca2SBarry Smith   x = (PetscReal)n;
1602b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1603d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1604b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1605d763cef2SBarry Smith   }
1606d763cef2SBarry Smith   PetscFunctionReturn(0);
1607d763cef2SBarry Smith }
1608d763cef2SBarry Smith 
16094a2ae208SSatish Balay #undef __FUNCT__
16104a2ae208SSatish Balay #define __FUNCT__ "TSLGMonitorDestroy"
1611d763cef2SBarry Smith /*@C
1612d763cef2SBarry Smith    TSLGMonitorDestroy - Destroys a line graph context that was created
1613d763cef2SBarry Smith    with TSLGMonitorCreate().
1614d763cef2SBarry Smith 
1615b0a32e0cSBarry Smith    Collective on PetscDrawLG
1616d763cef2SBarry Smith 
1617d763cef2SBarry Smith    Input Parameter:
1618d763cef2SBarry Smith .  draw - the drawing context
1619d763cef2SBarry Smith 
1620d763cef2SBarry Smith    Level: intermediate
1621d763cef2SBarry Smith 
1622d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1623d763cef2SBarry Smith 
1624d763cef2SBarry Smith .seealso: TSLGMonitorCreate(),  TSSetMonitor(), TSLGMonitor();
1625d763cef2SBarry Smith @*/
1626b0a32e0cSBarry Smith int TSLGMonitorDestroy(PetscDrawLG drawlg)
1627d763cef2SBarry Smith {
1628b0a32e0cSBarry Smith   PetscDraw draw;
1629d763cef2SBarry Smith   int       ierr;
1630d763cef2SBarry Smith 
1631d763cef2SBarry Smith   PetscFunctionBegin;
1632b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1633b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1634b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1635d763cef2SBarry Smith   PetscFunctionReturn(0);
1636d763cef2SBarry Smith }
1637d763cef2SBarry Smith 
16384a2ae208SSatish Balay #undef __FUNCT__
16394a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1640d763cef2SBarry Smith /*@
1641d763cef2SBarry Smith    TSGetTime - Gets the current time.
1642d763cef2SBarry Smith 
1643d763cef2SBarry Smith    Not Collective
1644d763cef2SBarry Smith 
1645d763cef2SBarry Smith    Input Parameter:
1646d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1647d763cef2SBarry Smith 
1648d763cef2SBarry Smith    Output Parameter:
1649d763cef2SBarry Smith .  t  - the current time
1650d763cef2SBarry Smith 
1651d763cef2SBarry Smith    Contributed by: Matthew Knepley
1652d763cef2SBarry Smith 
1653d763cef2SBarry Smith    Level: beginner
1654d763cef2SBarry Smith 
1655d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1656d763cef2SBarry Smith 
1657d763cef2SBarry Smith .keywords: TS, get, time
1658d763cef2SBarry Smith @*/
165987828ca2SBarry Smith int TSGetTime(TS ts,PetscReal* t)
1660d763cef2SBarry Smith {
1661d763cef2SBarry Smith   PetscFunctionBegin;
1662d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1663d763cef2SBarry Smith   *t = ts->ptime;
1664d763cef2SBarry Smith   PetscFunctionReturn(0);
1665d763cef2SBarry Smith }
1666d763cef2SBarry Smith 
16674a2ae208SSatish Balay #undef __FUNCT__
16684a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1669d763cef2SBarry Smith /*@C
1670d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1671d763cef2SBarry Smith    TS options in the database.
1672d763cef2SBarry Smith 
1673d763cef2SBarry Smith    Collective on TS
1674d763cef2SBarry Smith 
1675d763cef2SBarry Smith    Input Parameter:
1676d763cef2SBarry Smith +  ts     - The TS context
1677d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1678d763cef2SBarry Smith 
1679d763cef2SBarry Smith    Notes:
1680d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1681d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1682d763cef2SBarry Smith    hyphen.
1683d763cef2SBarry Smith 
1684d763cef2SBarry Smith    Contributed by: Matthew Knepley
1685d763cef2SBarry Smith 
1686d763cef2SBarry Smith    Level: advanced
1687d763cef2SBarry Smith 
1688d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1689d763cef2SBarry Smith 
1690d763cef2SBarry Smith .seealso: TSSetFromOptions()
1691d763cef2SBarry Smith 
1692d763cef2SBarry Smith @*/
1693d763cef2SBarry Smith int TSSetOptionsPrefix(TS ts,char *prefix)
1694d763cef2SBarry Smith {
1695d763cef2SBarry Smith   int ierr;
1696d763cef2SBarry Smith 
1697d763cef2SBarry Smith   PetscFunctionBegin;
1698d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1699d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1700d763cef2SBarry Smith   switch(ts->problem_type) {
1701d763cef2SBarry Smith     case TS_NONLINEAR:
1702d763cef2SBarry Smith       ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1703d763cef2SBarry Smith       break;
1704d763cef2SBarry Smith     case TS_LINEAR:
1705d763cef2SBarry Smith       ierr = SLESSetOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr);
1706d763cef2SBarry Smith       break;
1707d763cef2SBarry Smith   }
1708d763cef2SBarry Smith   PetscFunctionReturn(0);
1709d763cef2SBarry Smith }
1710d763cef2SBarry Smith 
1711d763cef2SBarry Smith 
17124a2ae208SSatish Balay #undef __FUNCT__
17134a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1714d763cef2SBarry Smith /*@C
1715d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1716d763cef2SBarry Smith    TS options in the database.
1717d763cef2SBarry Smith 
1718d763cef2SBarry Smith    Collective on TS
1719d763cef2SBarry Smith 
1720d763cef2SBarry Smith    Input Parameter:
1721d763cef2SBarry Smith +  ts     - The TS context
1722d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1723d763cef2SBarry Smith 
1724d763cef2SBarry Smith    Notes:
1725d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1726d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1727d763cef2SBarry Smith    hyphen.
1728d763cef2SBarry Smith 
1729d763cef2SBarry Smith    Contributed by: Matthew Knepley
1730d763cef2SBarry Smith 
1731d763cef2SBarry Smith    Level: advanced
1732d763cef2SBarry Smith 
1733d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1734d763cef2SBarry Smith 
1735d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1736d763cef2SBarry Smith 
1737d763cef2SBarry Smith @*/
1738d763cef2SBarry Smith int TSAppendOptionsPrefix(TS ts,char *prefix)
1739d763cef2SBarry Smith {
1740d763cef2SBarry Smith   int ierr;
1741d763cef2SBarry Smith 
1742d763cef2SBarry Smith   PetscFunctionBegin;
1743d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1744d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1745d763cef2SBarry Smith   switch(ts->problem_type) {
1746d763cef2SBarry Smith     case TS_NONLINEAR:
1747d763cef2SBarry Smith       ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
1748d763cef2SBarry Smith       break;
1749d763cef2SBarry Smith     case TS_LINEAR:
1750d763cef2SBarry Smith       ierr = SLESAppendOptionsPrefix(ts->sles,prefix);CHKERRQ(ierr);
1751d763cef2SBarry Smith       break;
1752d763cef2SBarry Smith   }
1753d763cef2SBarry Smith   PetscFunctionReturn(0);
1754d763cef2SBarry Smith }
1755d763cef2SBarry Smith 
17564a2ae208SSatish Balay #undef __FUNCT__
17574a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1758d763cef2SBarry Smith /*@C
1759d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1760d763cef2SBarry Smith    TS options in the database.
1761d763cef2SBarry Smith 
1762d763cef2SBarry Smith    Not Collective
1763d763cef2SBarry Smith 
1764d763cef2SBarry Smith    Input Parameter:
1765d763cef2SBarry Smith .  ts - The TS context
1766d763cef2SBarry Smith 
1767d763cef2SBarry Smith    Output Parameter:
1768d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1769d763cef2SBarry Smith 
1770d763cef2SBarry Smith    Contributed by: Matthew Knepley
1771d763cef2SBarry Smith 
1772d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1773d763cef2SBarry Smith    sufficient length to hold the prefix.
1774d763cef2SBarry Smith 
1775d763cef2SBarry Smith    Level: intermediate
1776d763cef2SBarry Smith 
1777d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1778d763cef2SBarry Smith 
1779d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1780d763cef2SBarry Smith @*/
1781d763cef2SBarry Smith int TSGetOptionsPrefix(TS ts,char **prefix)
1782d763cef2SBarry Smith {
1783d763cef2SBarry Smith   int ierr;
1784d763cef2SBarry Smith 
1785d763cef2SBarry Smith   PetscFunctionBegin;
1786d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1787d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1788d763cef2SBarry Smith   PetscFunctionReturn(0);
1789d763cef2SBarry Smith }
1790d763cef2SBarry Smith 
17914a2ae208SSatish Balay #undef __FUNCT__
17924a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSMatrix"
1793d763cef2SBarry Smith /*@C
1794d763cef2SBarry Smith    TSGetRHSMatrix - Returns the matrix A at the present timestep.
1795d763cef2SBarry Smith 
1796d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1797d763cef2SBarry Smith 
1798d763cef2SBarry Smith    Input Parameter:
1799d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1800d763cef2SBarry Smith 
1801d763cef2SBarry Smith    Output Parameters:
1802d763cef2SBarry Smith +  A   - The matrix A, where U_t = A(t) U
1803d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as A
1804d763cef2SBarry Smith -  ctx - User-defined context for matrix evaluation routine
1805d763cef2SBarry Smith 
1806d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1807d763cef2SBarry Smith 
1808d763cef2SBarry Smith    Contributed by: Matthew Knepley
1809d763cef2SBarry Smith 
1810d763cef2SBarry Smith    Level: intermediate
1811d763cef2SBarry Smith 
1812d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
1813d763cef2SBarry Smith 
1814d763cef2SBarry Smith .keywords: TS, timestep, get, matrix
1815d763cef2SBarry Smith 
1816d763cef2SBarry Smith @*/
1817d763cef2SBarry Smith int TSGetRHSMatrix(TS ts,Mat *A,Mat *M,void **ctx)
1818d763cef2SBarry Smith {
1819d763cef2SBarry Smith   PetscFunctionBegin;
1820d763cef2SBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE);
1821d763cef2SBarry Smith   if (A)   *A = ts->A;
1822d763cef2SBarry Smith   if (M)   *M = ts->B;
1823d763cef2SBarry Smith   if (ctx) *ctx = ts->jacP;
1824d763cef2SBarry Smith   PetscFunctionReturn(0);
1825d763cef2SBarry Smith }
1826d763cef2SBarry Smith 
18274a2ae208SSatish Balay #undef __FUNCT__
18284a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1829d763cef2SBarry Smith /*@C
1830d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1831d763cef2SBarry Smith 
1832d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1833d763cef2SBarry Smith 
1834d763cef2SBarry Smith    Input Parameter:
1835d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1836d763cef2SBarry Smith 
1837d763cef2SBarry Smith    Output Parameters:
1838d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1839d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1840d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1841d763cef2SBarry Smith 
1842d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1843d763cef2SBarry Smith 
1844d763cef2SBarry Smith    Contributed by: Matthew Knepley
1845d763cef2SBarry Smith 
1846d763cef2SBarry Smith    Level: intermediate
1847d763cef2SBarry Smith 
1848d763cef2SBarry Smith .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber()
1849d763cef2SBarry Smith 
1850d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1851d763cef2SBarry Smith @*/
1852d763cef2SBarry Smith int TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1853d763cef2SBarry Smith {
1854d763cef2SBarry Smith   int ierr;
1855d763cef2SBarry Smith 
1856d763cef2SBarry Smith   PetscFunctionBegin;
1857d763cef2SBarry Smith   ierr = TSGetRHSMatrix(ts,J,M,ctx);CHKERRQ(ierr);
1858d763cef2SBarry Smith   PetscFunctionReturn(0);
1859d763cef2SBarry Smith }
1860d763cef2SBarry Smith 
1861d763cef2SBarry Smith /*MC
1862f1af5d2fSBarry Smith    TSRegisterDynamic - Adds a method to the timestepping solver package.
1863d763cef2SBarry Smith 
1864d763cef2SBarry Smith    Synopsis:
18653a7fca6bSBarry Smith    int TSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(TS))
1866d763cef2SBarry Smith 
1867d763cef2SBarry Smith    Not collective
1868d763cef2SBarry Smith 
1869d763cef2SBarry Smith    Input Parameters:
1870d763cef2SBarry Smith +  name_solver - name of a new user-defined solver
1871d763cef2SBarry Smith .  path - path (either absolute or relative) the library containing this solver
1872d763cef2SBarry Smith .  name_create - name of routine to create method context
1873d763cef2SBarry Smith -  routine_create - routine to create method context
1874d763cef2SBarry Smith 
1875d763cef2SBarry Smith    Notes:
1876f1af5d2fSBarry Smith    TSRegisterDynamic() may be called multiple times to add several user-defined solvers.
1877d763cef2SBarry Smith 
1878d763cef2SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1879d763cef2SBarry Smith    is ignored.
1880d763cef2SBarry Smith 
1881d763cef2SBarry Smith    Sample usage:
1882d763cef2SBarry Smith .vb
1883f1af5d2fSBarry Smith    TSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
1884d763cef2SBarry Smith               "MySolverCreate",MySolverCreate);
1885d763cef2SBarry Smith .ve
1886d763cef2SBarry Smith 
1887d763cef2SBarry Smith    Then, your solver can be chosen with the procedural interface via
1888d763cef2SBarry Smith $     TSSetType(ts,"my_solver")
1889d763cef2SBarry Smith    or at runtime via the option
1890d763cef2SBarry Smith $     -ts_type my_solver
1891d763cef2SBarry Smith 
1892d763cef2SBarry Smith    Level: advanced
1893d763cef2SBarry Smith 
18942aeb12f1SSatish Balay    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
18955ba88a07SLois Curfman McInnes    and others of the form ${any_environmental_variable} occuring in pathname will be
18965ba88a07SLois Curfman McInnes    replaced with appropriate values.
1897d763cef2SBarry Smith 
1898d763cef2SBarry Smith .keywords: TS, register
1899d763cef2SBarry Smith 
1900d763cef2SBarry Smith .seealso: TSRegisterAll(), TSRegisterDestroy()
1901d763cef2SBarry Smith M*/
1902d763cef2SBarry Smith 
19034a2ae208SSatish Balay #undef __FUNCT__
19044a2ae208SSatish Balay #define __FUNCT__ "TSRegister"
1905000e7ae3SMatthew Knepley int TSRegister(const char sname[], const char path[], const char name[], int (*function)(TS))
1906d763cef2SBarry Smith {
1907d763cef2SBarry Smith   char fullname[256];
1908549d3d68SSatish Balay   int  ierr;
1909d763cef2SBarry Smith 
1910d763cef2SBarry Smith   PetscFunctionBegin;
1911b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1912c134de8dSSatish Balay   ierr = PetscFListAdd(&TSList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1913d763cef2SBarry Smith   PetscFunctionReturn(0);
1914d763cef2SBarry Smith }
19151713a123SBarry Smith 
19161713a123SBarry Smith #undef __FUNCT__
19171713a123SBarry Smith #define __FUNCT__ "TSVecViewMonitor"
19181713a123SBarry Smith /*@C
19191713a123SBarry Smith    TSVecViewMonitor - Monitors progress of the TS solvers by calling
19201713a123SBarry Smith    VecView() for the solution at each timestep
19211713a123SBarry Smith 
19221713a123SBarry Smith    Collective on TS
19231713a123SBarry Smith 
19241713a123SBarry Smith    Input Parameters:
19251713a123SBarry Smith +  ts - the TS context
19261713a123SBarry Smith .  step - current time-step
1927142b95e3SSatish Balay .  ptime - current time
19281713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
19291713a123SBarry Smith 
19301713a123SBarry Smith    Level: intermediate
19311713a123SBarry Smith 
19321713a123SBarry Smith .keywords: TS,  vector, monitor, view
19331713a123SBarry Smith 
19341713a123SBarry Smith .seealso: TSSetMonitor(), TSDefaultMonitor(), VecView()
19351713a123SBarry Smith @*/
1936142b95e3SSatish Balay int TSVecViewMonitor(TS ts,int step,PetscReal ptime,Vec x,void *dummy)
19371713a123SBarry Smith {
19381713a123SBarry Smith   int         ierr;
19391713a123SBarry Smith   PetscViewer viewer = (PetscViewer) dummy;
19401713a123SBarry Smith 
19411713a123SBarry Smith   PetscFunctionBegin;
19421713a123SBarry Smith   if (!viewer) {
19431713a123SBarry Smith     MPI_Comm comm;
19441713a123SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
19451713a123SBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
19461713a123SBarry Smith   }
19471713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
19481713a123SBarry Smith   PetscFunctionReturn(0);
19491713a123SBarry Smith }
19501713a123SBarry Smith 
19511713a123SBarry Smith 
19521713a123SBarry Smith 
1953