xref: /petsc/src/ts/interface/ts.c (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
163dd3a1aSKris Buschelman #define PETSCTS_DLL
263dd3a1aSKris Buschelman 
37c4f633dSBarry Smith #include "private/tsimpl.h"        /*I "petscts.h"  I*/
4d763cef2SBarry Smith 
5d5ba7fb7SMatthew Knepley /* Logging support */
6*0700a824SBarry Smith PetscClassId PETSCTS_DLLEXPORT TS_CLASSID;
7166c7f25SBarry Smith PetscLogEvent  TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
8d405a339SMatthew Knepley 
94a2ae208SSatish Balay #undef __FUNCT__
10bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
11bdad233fSMatthew Knepley /*
12bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
13bdad233fSMatthew Knepley 
14bdad233fSMatthew Knepley   Collective on TS
15bdad233fSMatthew Knepley 
16bdad233fSMatthew Knepley   Input Parameter:
17bdad233fSMatthew Knepley . ts - The ts
18bdad233fSMatthew Knepley 
19bdad233fSMatthew Knepley   Level: intermediate
20bdad233fSMatthew Knepley 
21bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
22bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
23bdad233fSMatthew Knepley */
246849ba73SBarry Smith static PetscErrorCode TSSetTypeFromOptions(TS ts)
25bdad233fSMatthew Knepley {
26bdad233fSMatthew Knepley   PetscTruth     opt;
272fc52814SBarry Smith   const char     *defaultType;
28bdad233fSMatthew Knepley   char           typeName[256];
29dfbe8321SBarry Smith   PetscErrorCode ierr;
30bdad233fSMatthew Knepley 
31bdad233fSMatthew Knepley   PetscFunctionBegin;
327adad957SLisandro Dalcin   if (((PetscObject)ts)->type_name) {
337adad957SLisandro Dalcin     defaultType = ((PetscObject)ts)->type_name;
34bdad233fSMatthew Knepley   } else {
359596e0b4SJed Brown     defaultType = TSEULER;
36bdad233fSMatthew Knepley   }
37bdad233fSMatthew Knepley 
38cce0b1b2SLisandro Dalcin   if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
39bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
40a7cc72afSBarry Smith   if (opt) {
41bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
42bdad233fSMatthew Knepley   } else {
43bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
44bdad233fSMatthew Knepley   }
45bdad233fSMatthew Knepley   PetscFunctionReturn(0);
46bdad233fSMatthew Knepley }
47bdad233fSMatthew Knepley 
48bdad233fSMatthew Knepley #undef __FUNCT__
49bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
50bdad233fSMatthew Knepley /*@
51bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
52bdad233fSMatthew Knepley 
53bdad233fSMatthew Knepley    Collective on TS
54bdad233fSMatthew Knepley 
55bdad233fSMatthew Knepley    Input Parameter:
56bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
57bdad233fSMatthew Knepley 
58bdad233fSMatthew Knepley    Options Database Keys:
599596e0b4SJed Brown +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCRANK_NICHOLSON
60bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
61bdad233fSMatthew Knepley .  -ts_max_time time - maximum time to compute to
62bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
63bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
64a6570f20SBarry Smith -  -ts_monitor_draw - plot information at each timestep
65bdad233fSMatthew Knepley 
66bdad233fSMatthew Knepley    Level: beginner
67bdad233fSMatthew Knepley 
68bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
69bdad233fSMatthew Knepley 
70a313700dSBarry Smith .seealso: TSGetType()
71bdad233fSMatthew Knepley @*/
7263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetFromOptions(TS ts)
73bdad233fSMatthew Knepley {
74bdad233fSMatthew Knepley   PetscReal               dt;
75eabae89aSBarry Smith   PetscTruth              opt,flg;
76dfbe8321SBarry Smith   PetscErrorCode          ierr;
77a34d58ebSBarry Smith   PetscViewerASCIIMonitor monviewer;
78eabae89aSBarry Smith   char                    monfilename[PETSC_MAX_PATH_LEN];
79bdad233fSMatthew Knepley 
80bdad233fSMatthew Knepley   PetscFunctionBegin;
81*0700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
827adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr);
83bdad233fSMatthew Knepley 
84bdad233fSMatthew Knepley     /* Handle generic TS options */
85bdad233fSMatthew Knepley     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
86bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
87bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
88bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
89a7cc72afSBarry Smith     if (opt) {
90bdad233fSMatthew Knepley       ts->initial_time_step = ts->time_step = dt;
91bdad233fSMatthew Knepley     }
92bdad233fSMatthew Knepley 
93bdad233fSMatthew Knepley     /* Monitor options */
94a6570f20SBarry Smith     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
95eabae89aSBarry Smith     if (flg) {
96050a712dSBarry Smith       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,monfilename,((PetscObject)ts)->tablevel,&monviewer);CHKERRQ(ierr);
97a34d58ebSBarry Smith       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
98bdad233fSMatthew Knepley     }
9990d69ab7SBarry Smith     opt  = PETSC_FALSE;
10090d69ab7SBarry Smith     ierr = PetscOptionsTruth("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
101a7cc72afSBarry Smith     if (opt) {
102a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
103bdad233fSMatthew Knepley     }
10490d69ab7SBarry Smith     opt  = PETSC_FALSE;
10590d69ab7SBarry Smith     ierr = PetscOptionsTruth("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
106a7cc72afSBarry Smith     if (opt) {
107a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
108bdad233fSMatthew Knepley     }
109bdad233fSMatthew Knepley 
110bdad233fSMatthew Knepley     /* Handle TS type options */
111bdad233fSMatthew Knepley     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
112bdad233fSMatthew Knepley 
113bdad233fSMatthew Knepley     /* Handle specific TS options */
114abc0a331SBarry Smith     if (ts->ops->setfromoptions) {
115bdad233fSMatthew Knepley       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
116bdad233fSMatthew Knepley     }
117bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
118bdad233fSMatthew Knepley 
119bdad233fSMatthew Knepley   /* Handle subobject options */
120bdad233fSMatthew Knepley   switch(ts->problem_type) {
121156fc9a6SMatthew Knepley     /* Should check for implicit/explicit */
122bdad233fSMatthew Knepley   case TS_LINEAR:
123abc0a331SBarry Smith     if (ts->ksp) {
1248beb423aSHong Zhang       ierr = KSPSetOperators(ts->ksp,ts->Arhs,ts->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
12594b7f48cSBarry Smith       ierr = KSPSetFromOptions(ts->ksp);CHKERRQ(ierr);
126156fc9a6SMatthew Knepley     }
127bdad233fSMatthew Knepley     break;
128bdad233fSMatthew Knepley   case TS_NONLINEAR:
129abc0a331SBarry Smith     if (ts->snes) {
1307c236d22SBarry Smith       /* this is a bit of a hack, but it gets the matrix information into SNES earlier
1317c236d22SBarry Smith          so that SNES and KSP have more information to pick reasonable defaults
1327c0b301bSJed Brown          before they allow users to set options
1337c0b301bSJed Brown        * If ts->A has been set at this point, we are probably using the implicit form
1347c0b301bSJed Brown          and Arhs will never be used. */
1357c0b301bSJed Brown       ierr = SNESSetJacobian(ts->snes,ts->A?ts->A:ts->Arhs,ts->B,0,ts);CHKERRQ(ierr);
136bdad233fSMatthew Knepley       ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
137156fc9a6SMatthew Knepley     }
138bdad233fSMatthew Knepley     break;
139bdad233fSMatthew Knepley   default:
14077431f27SBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type);
141bdad233fSMatthew Knepley   }
142bdad233fSMatthew Knepley 
143bdad233fSMatthew Knepley   PetscFunctionReturn(0);
144bdad233fSMatthew Knepley }
145bdad233fSMatthew Knepley 
146bdad233fSMatthew Knepley #undef  __FUNCT__
147bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
148bdad233fSMatthew Knepley /*@
149bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
150bdad233fSMatthew Knepley 
151bdad233fSMatthew Knepley   Collective on TS
152bdad233fSMatthew Knepley 
153bdad233fSMatthew Knepley   Input Parameter:
154bdad233fSMatthew Knepley . ts - The ts
155bdad233fSMatthew Knepley 
156bdad233fSMatthew Knepley   Level: intermediate
157bdad233fSMatthew Knepley 
158bdad233fSMatthew Knepley .keywords: TS, view, options, database
159bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
160bdad233fSMatthew Knepley @*/
16163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSViewFromOptions(TS ts,const char title[])
162bdad233fSMatthew Knepley {
163bdad233fSMatthew Knepley   PetscViewer    viewer;
164bdad233fSMatthew Knepley   PetscDraw      draw;
16590d69ab7SBarry Smith   PetscTruth     opt = PETSC_FALSE;
166e10c49a3SBarry Smith   char           fileName[PETSC_MAX_PATH_LEN];
167dfbe8321SBarry Smith   PetscErrorCode ierr;
168bdad233fSMatthew Knepley 
169bdad233fSMatthew Knepley   PetscFunctionBegin;
1707adad957SLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr);
171eabae89aSBarry Smith   if (opt && !PetscPreLoadingOn) {
1727adad957SLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr);
173bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
174bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
175bdad233fSMatthew Knepley   }
1768e83347fSKai Germaschewski   opt = PETSC_FALSE;
17790d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr);
178a7cc72afSBarry Smith   if (opt) {
1797adad957SLisandro Dalcin     ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
180bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
181a7cc72afSBarry Smith     if (title) {
1821836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
183bdad233fSMatthew Knepley     } else {
184bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr);
1857adad957SLisandro Dalcin       ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr);
186bdad233fSMatthew Knepley     }
187bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
188bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
189bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
190bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
191bdad233fSMatthew Knepley   }
192bdad233fSMatthew Knepley   PetscFunctionReturn(0);
193bdad233fSMatthew Knepley }
194bdad233fSMatthew Knepley 
195bdad233fSMatthew Knepley #undef __FUNCT__
1964a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
197a7a1495cSBarry Smith /*@
1988c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
199a7a1495cSBarry Smith       set with TSSetRHSJacobian().
200a7a1495cSBarry Smith 
201a7a1495cSBarry Smith    Collective on TS and Vec
202a7a1495cSBarry Smith 
203a7a1495cSBarry Smith    Input Parameters:
204316643e7SJed Brown +  ts - the TS context
205a7a1495cSBarry Smith .  t - current timestep
206a7a1495cSBarry Smith -  x - input vector
207a7a1495cSBarry Smith 
208a7a1495cSBarry Smith    Output Parameters:
209a7a1495cSBarry Smith +  A - Jacobian matrix
210a7a1495cSBarry Smith .  B - optional preconditioning matrix
211a7a1495cSBarry Smith -  flag - flag indicating matrix structure
212a7a1495cSBarry Smith 
213a7a1495cSBarry Smith    Notes:
214a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
215a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
216a7a1495cSBarry Smith 
21794b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
218a7a1495cSBarry Smith    flag parameter.
219a7a1495cSBarry Smith 
220a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
221a7a1495cSBarry Smith 
222a7a1495cSBarry Smith    Level: developer
223a7a1495cSBarry Smith 
224a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
225a7a1495cSBarry Smith 
22694b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
227a7a1495cSBarry Smith @*/
22863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
229a7a1495cSBarry Smith {
230dfbe8321SBarry Smith   PetscErrorCode ierr;
231a7a1495cSBarry Smith 
232a7a1495cSBarry Smith   PetscFunctionBegin;
233*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
234*0700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
235c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
236a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
23729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
238a7a1495cSBarry Smith   }
239000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
240d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
241a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
242a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
243000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
244a7a1495cSBarry Smith     PetscStackPop;
245d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
246a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
247*0700a824SBarry Smith     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
248*0700a824SBarry Smith     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
249ef66eb69SBarry Smith   } else {
250ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
251ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252ef66eb69SBarry Smith     if (*A != *B) {
253ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
254ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255ef66eb69SBarry Smith     }
256ef66eb69SBarry Smith   }
257a7a1495cSBarry Smith   PetscFunctionReturn(0);
258a7a1495cSBarry Smith }
259a7a1495cSBarry Smith 
2604a2ae208SSatish Balay #undef __FUNCT__
2614a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
262316643e7SJed Brown /*@
263d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
264d763cef2SBarry Smith 
265316643e7SJed Brown    Collective on TS and Vec
266316643e7SJed Brown 
267316643e7SJed Brown    Input Parameters:
268316643e7SJed Brown +  ts - the TS context
269316643e7SJed Brown .  t - current time
270316643e7SJed Brown -  x - state vector
271316643e7SJed Brown 
272316643e7SJed Brown    Output Parameter:
273316643e7SJed Brown .  y - right hand side
274316643e7SJed Brown 
275316643e7SJed Brown    Note:
276316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
277316643e7SJed Brown    is used internally within the nonlinear solvers.
278316643e7SJed Brown 
279316643e7SJed Brown    If the user did not provide a function but merely a matrix,
280d763cef2SBarry Smith    this routine applies the matrix.
281316643e7SJed Brown 
282316643e7SJed Brown    Level: developer
283316643e7SJed Brown 
284316643e7SJed Brown .keywords: TS, compute
285316643e7SJed Brown 
286316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction()
287316643e7SJed Brown @*/
288dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
289d763cef2SBarry Smith {
290dfbe8321SBarry Smith   PetscErrorCode ierr;
291d763cef2SBarry Smith 
292d763cef2SBarry Smith   PetscFunctionBegin;
293*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
294*0700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
295*0700a824SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
296d763cef2SBarry Smith 
297d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
298000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
299d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
300000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
301d763cef2SBarry Smith     PetscStackPop;
3027c922b88SBarry Smith   } else {
303000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
304d763cef2SBarry Smith       MatStructure flg;
305d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
3068beb423aSHong Zhang       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
307d763cef2SBarry Smith       PetscStackPop;
308d763cef2SBarry Smith     }
3098beb423aSHong Zhang     ierr = MatMult(ts->Arhs,x,y);CHKERRQ(ierr);
3107c922b88SBarry Smith   }
311d763cef2SBarry Smith 
312d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
313d763cef2SBarry Smith 
314d763cef2SBarry Smith   PetscFunctionReturn(0);
315d763cef2SBarry Smith }
316d763cef2SBarry Smith 
3174a2ae208SSatish Balay #undef __FUNCT__
318316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction"
319316643e7SJed Brown /*@
320316643e7SJed Brown    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
321316643e7SJed Brown 
322316643e7SJed Brown    Collective on TS and Vec
323316643e7SJed Brown 
324316643e7SJed Brown    Input Parameters:
325316643e7SJed Brown +  ts - the TS context
326316643e7SJed Brown .  t - current time
327316643e7SJed Brown .  X - state vector
328316643e7SJed Brown -  Xdot - time derivative of state vector
329316643e7SJed Brown 
330316643e7SJed Brown    Output Parameter:
331316643e7SJed Brown .  Y - right hand side
332316643e7SJed Brown 
333316643e7SJed Brown    Note:
334316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
335316643e7SJed Brown    is used internally within the nonlinear solvers.
336316643e7SJed Brown 
337316643e7SJed Brown    If the user did did not write their equations in implicit form, this
338316643e7SJed Brown    function recasts them in implicit form.
339316643e7SJed Brown 
340316643e7SJed Brown    Level: developer
341316643e7SJed Brown 
342316643e7SJed Brown .keywords: TS, compute
343316643e7SJed Brown 
344316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction()
345316643e7SJed Brown @*/
346316643e7SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y)
347316643e7SJed Brown {
348316643e7SJed Brown   PetscErrorCode ierr;
349316643e7SJed Brown 
350316643e7SJed Brown   PetscFunctionBegin;
351*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
352*0700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
353*0700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
354*0700a824SBarry Smith   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
355316643e7SJed Brown 
356316643e7SJed Brown   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
357316643e7SJed Brown   if (ts->ops->ifunction) {
358316643e7SJed Brown     PetscStackPush("TS user implicit function");
359316643e7SJed Brown     ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr);
360316643e7SJed Brown     PetscStackPop;
361316643e7SJed Brown   } else {
362316643e7SJed Brown     if (ts->ops->rhsfunction) {
363316643e7SJed Brown       PetscStackPush("TS user right-hand-side function");
364316643e7SJed Brown       ierr = (*ts->ops->rhsfunction)(ts,t,X,Y,ts->funP);CHKERRQ(ierr);
365316643e7SJed Brown       PetscStackPop;
366316643e7SJed Brown     } else {
367316643e7SJed Brown       if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
368316643e7SJed Brown         MatStructure flg;
3694a6899ffSJed Brown         /* Note: flg is not being used.
3704a6899ffSJed Brown            For it to be useful, we'd have to cache it and then apply it in TSComputeIJacobian.
3714a6899ffSJed Brown         */
372316643e7SJed Brown         PetscStackPush("TS user right-hand-side matrix function");
373316643e7SJed Brown         ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
374316643e7SJed Brown         PetscStackPop;
375316643e7SJed Brown       }
376316643e7SJed Brown       ierr = MatMult(ts->Arhs,X,Y);CHKERRQ(ierr);
377316643e7SJed Brown     }
378316643e7SJed Brown 
3794a6899ffSJed Brown     /* Convert to implicit form: F(X,Xdot) = Alhs * Xdot - Frhs(X) */
3804a6899ffSJed Brown     if (ts->Alhs) {
3814a6899ffSJed Brown       if (ts->ops->lhsmatrix) {
3824a6899ffSJed Brown         MatStructure flg;
3834a6899ffSJed Brown         PetscStackPush("TS user left-hand-side matrix function");
3844a6899ffSJed Brown         ierr = (*ts->ops->lhsmatrix)(ts,t,&ts->Alhs,PETSC_NULL,&flg,ts->jacP);CHKERRQ(ierr);
3854a6899ffSJed Brown         PetscStackPop;
3864a6899ffSJed Brown       }
3874a6899ffSJed Brown       ierr = VecScale(Y,-1.);CHKERRQ(ierr);
3884a6899ffSJed Brown       ierr = MatMultAdd(ts->Alhs,Xdot,Y,Y);CHKERRQ(ierr);
3894a6899ffSJed Brown     } else {
390ace68cafSJed Brown       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
391316643e7SJed Brown     }
3924a6899ffSJed Brown   }
393316643e7SJed Brown   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
394316643e7SJed Brown   PetscFunctionReturn(0);
395316643e7SJed Brown }
396316643e7SJed Brown 
397316643e7SJed Brown #undef __FUNCT__
398316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian"
399316643e7SJed Brown /*@
400316643e7SJed Brown    TSComputeIJacobian - Evaluates the Jacobian of the DAE
401316643e7SJed Brown 
402316643e7SJed Brown    Collective on TS and Vec
403316643e7SJed Brown 
404316643e7SJed Brown    Input
405316643e7SJed Brown       Input Parameters:
406316643e7SJed Brown +  ts - the TS context
407316643e7SJed Brown .  t - current timestep
408316643e7SJed Brown .  X - state vector
409316643e7SJed Brown .  Xdot - time derivative of state vector
410316643e7SJed Brown -  shift - shift to apply, see note below
411316643e7SJed Brown 
412316643e7SJed Brown    Output Parameters:
413316643e7SJed Brown +  A - Jacobian matrix
414316643e7SJed Brown .  B - optional preconditioning matrix
415316643e7SJed Brown -  flag - flag indicating matrix structure
416316643e7SJed Brown 
417316643e7SJed Brown    Notes:
418316643e7SJed Brown    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
419316643e7SJed Brown 
420316643e7SJed Brown    dF/dX + shift*dF/dXdot
421316643e7SJed Brown 
422316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
423316643e7SJed Brown    is used internally within the nonlinear solvers.
424316643e7SJed Brown 
425316643e7SJed Brown    TSComputeIJacobian() is valid only for TS_NONLINEAR
426316643e7SJed Brown 
427316643e7SJed Brown    Level: developer
428316643e7SJed Brown 
429316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix
430316643e7SJed Brown 
431316643e7SJed Brown .seealso:  TSSetIJacobian()
43263495f91SJed Brown @*/
433316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg)
434316643e7SJed Brown {
435316643e7SJed Brown   PetscErrorCode ierr;
436316643e7SJed Brown 
437316643e7SJed Brown   PetscFunctionBegin;
438*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
439*0700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
440*0700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
441316643e7SJed Brown   PetscValidPointer(A,6);
442*0700a824SBarry Smith   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
443316643e7SJed Brown   PetscValidPointer(B,7);
444*0700a824SBarry Smith   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
445316643e7SJed Brown   PetscValidPointer(flg,8);
446316643e7SJed Brown 
4474a6899ffSJed Brown   *flg = SAME_NONZERO_PATTERN;  /* In case it we're solving a linear problem in which case it wouldn't get initialized below. */
448316643e7SJed Brown   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
449316643e7SJed Brown   if (ts->ops->ijacobian) {
450316643e7SJed Brown     PetscStackPush("TS user implicit Jacobian");
451316643e7SJed Brown     ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr);
452316643e7SJed Brown     PetscStackPop;
453316643e7SJed Brown   } else {
454316643e7SJed Brown     if (ts->ops->rhsjacobian) {
455316643e7SJed Brown       PetscStackPush("TS user right-hand-side Jacobian");
456316643e7SJed Brown       ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
457316643e7SJed Brown       PetscStackPop;
458316643e7SJed Brown     } else {
459316643e7SJed Brown       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
460316643e7SJed Brown       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
461316643e7SJed Brown       if (*A != *B) {
462316643e7SJed Brown         ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
463316643e7SJed Brown         ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
464316643e7SJed Brown       }
465316643e7SJed Brown     }
466316643e7SJed Brown 
467316643e7SJed Brown     /* Convert to implicit form */
468316643e7SJed Brown     /* inefficient because these operations will normally traverse all matrix elements twice */
469316643e7SJed Brown     ierr = MatScale(*A,-1);CHKERRQ(ierr);
4704a6899ffSJed Brown     if (ts->Alhs) {
4714a6899ffSJed Brown       ierr = MatAXPY(*A,shift,ts->Alhs,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
4724a6899ffSJed Brown     } else {
473316643e7SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
4744a6899ffSJed Brown     }
475316643e7SJed Brown     if (*A != *B) {
476316643e7SJed Brown       ierr = MatScale(*B,-1);CHKERRQ(ierr);
477316643e7SJed Brown       ierr = MatShift(*B,shift);CHKERRQ(ierr);
478316643e7SJed Brown     }
479316643e7SJed Brown   }
480316643e7SJed Brown   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
481316643e7SJed Brown   PetscFunctionReturn(0);
482316643e7SJed Brown }
483316643e7SJed Brown 
484316643e7SJed Brown #undef __FUNCT__
4854a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
486d763cef2SBarry Smith /*@C
487d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
488d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
489d763cef2SBarry Smith 
490d763cef2SBarry Smith     Collective on TS
491d763cef2SBarry Smith 
492d763cef2SBarry Smith     Input Parameters:
493d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
494d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
495d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
496d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
497d763cef2SBarry Smith 
498d763cef2SBarry Smith     Calling sequence of func:
49987828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
500d763cef2SBarry Smith 
501d763cef2SBarry Smith +   t - current timestep
502d763cef2SBarry Smith .   u - input vector
503d763cef2SBarry Smith .   F - function vector
504d763cef2SBarry Smith -   ctx - [optional] user-defined function context
505d763cef2SBarry Smith 
506d763cef2SBarry Smith     Important:
50776f2fa84SHong Zhang     The user MUST call either this routine or TSSetMatrices().
508d763cef2SBarry Smith 
509d763cef2SBarry Smith     Level: beginner
510d763cef2SBarry Smith 
511d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
512d763cef2SBarry Smith 
51376f2fa84SHong Zhang .seealso: TSSetMatrices()
514d763cef2SBarry Smith @*/
51563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
516d763cef2SBarry Smith {
517d763cef2SBarry Smith   PetscFunctionBegin;
518d763cef2SBarry Smith 
519*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
520d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
52129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
522d763cef2SBarry Smith   }
523000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
524d763cef2SBarry Smith   ts->funP             = ctx;
525d763cef2SBarry Smith   PetscFunctionReturn(0);
526d763cef2SBarry Smith }
527d763cef2SBarry Smith 
5284a2ae208SSatish Balay #undef __FUNCT__
52995f0b562SHong Zhang #define __FUNCT__ "TSSetMatrices"
53095f0b562SHong Zhang /*@C
53195f0b562SHong Zhang    TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs,
53295f0b562SHong Zhang    where Alhs(t) U_t = Arhs(t) U.
53395f0b562SHong Zhang 
53495f0b562SHong Zhang    Collective on TS
53595f0b562SHong Zhang 
53695f0b562SHong Zhang    Input Parameters:
53795f0b562SHong Zhang +  ts   - the TS context obtained from TSCreate()
53895f0b562SHong Zhang .  Arhs - matrix
53995f0b562SHong Zhang .  frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
54095f0b562SHong Zhang           if Arhs is not a function of t.
54195f0b562SHong Zhang .  Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix.
54295f0b562SHong Zhang .  flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
54395f0b562SHong Zhang           if Alhs is not a function of t.
54495f0b562SHong Zhang .  flag - flag indicating information about the matrix structure of Arhs and Alhs.
54595f0b562SHong Zhang           The available options are
54695f0b562SHong Zhang             SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs
54795f0b562SHong Zhang             DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs
54895f0b562SHong Zhang -  ctx  - [optional] user-defined context for private data for the
54995f0b562SHong Zhang           matrix evaluation routine (may be PETSC_NULL)
55095f0b562SHong Zhang 
55195f0b562SHong Zhang    Calling sequence of func:
55295f0b562SHong Zhang $     func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
55395f0b562SHong Zhang 
55495f0b562SHong Zhang +  t - current timestep
55595f0b562SHong Zhang .  A - matrix A, where U_t = A(t) U
55695f0b562SHong Zhang .  B - preconditioner matrix, usually the same as A
55795f0b562SHong Zhang .  flag - flag indicating information about the preconditioner matrix
55895f0b562SHong Zhang           structure (same as flag in KSPSetOperators())
55995f0b562SHong Zhang -  ctx - [optional] user-defined context for matrix evaluation routine
56095f0b562SHong Zhang 
56195f0b562SHong Zhang    Notes:
56295f0b562SHong Zhang    The routine func() takes Mat* as the matrix arguments rather than Mat.
56395f0b562SHong Zhang    This allows the matrix evaluation routine to replace Arhs or Alhs with a
56495f0b562SHong Zhang    completely new new matrix structure (not just different matrix elements)
56595f0b562SHong Zhang    when appropriate, for instance, if the nonzero structure is changing
56695f0b562SHong Zhang    throughout the global iterations.
56795f0b562SHong Zhang 
56895f0b562SHong Zhang    Important:
56995f0b562SHong Zhang    The user MUST call either this routine or TSSetRHSFunction().
57095f0b562SHong Zhang 
57195f0b562SHong Zhang    Level: beginner
57295f0b562SHong Zhang 
57395f0b562SHong Zhang .keywords: TS, timestep, set, matrix
57495f0b562SHong Zhang 
57595f0b562SHong Zhang .seealso: TSSetRHSFunction()
57695f0b562SHong Zhang @*/
57795f0b562SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSetMatrices(TS ts,Mat Arhs,PetscErrorCode (*frhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),Mat Alhs,PetscErrorCode (*flhs)(TS,PetscReal,Mat*,Mat*,MatStructure*,void*),MatStructure flag,void *ctx)
57895f0b562SHong Zhang {
57995f0b562SHong Zhang   PetscFunctionBegin;
580*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
58192af4f6aSHong Zhang   if (Arhs){
582*0700a824SBarry Smith     PetscValidHeaderSpecific(Arhs,MAT_CLASSID,2);
58395f0b562SHong Zhang     PetscCheckSameComm(ts,1,Arhs,2);
58495f0b562SHong Zhang     ts->Arhs           = Arhs;
58592af4f6aSHong Zhang     ts->ops->rhsmatrix = frhs;
58692af4f6aSHong Zhang   }
58792af4f6aSHong Zhang   if (Alhs){
588*0700a824SBarry Smith     PetscValidHeaderSpecific(Alhs,MAT_CLASSID,4);
58992af4f6aSHong Zhang     PetscCheckSameComm(ts,1,Alhs,4);
59095f0b562SHong Zhang     ts->Alhs           = Alhs;
59192af4f6aSHong Zhang     ts->ops->lhsmatrix = flhs;
59292af4f6aSHong Zhang   }
59392af4f6aSHong Zhang 
59492af4f6aSHong Zhang   ts->jacP           = ctx;
59595f0b562SHong Zhang   ts->matflg         = flag;
59695f0b562SHong Zhang   PetscFunctionReturn(0);
59795f0b562SHong Zhang }
598d763cef2SBarry Smith 
599aa644b49SHong Zhang #undef __FUNCT__
600cda39b92SHong Zhang #define __FUNCT__ "TSGetMatrices"
601cda39b92SHong Zhang /*@C
602cda39b92SHong Zhang    TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep,
603cda39b92SHong Zhang    where Alhs(t) U_t = Arhs(t) U.
604cda39b92SHong Zhang 
605cda39b92SHong Zhang    Not Collective, but parallel objects are returned if TS is parallel
606cda39b92SHong Zhang 
607cda39b92SHong Zhang    Input Parameter:
608cda39b92SHong Zhang .  ts  - The TS context obtained from TSCreate()
609cda39b92SHong Zhang 
610cda39b92SHong Zhang    Output Parameters:
611cda39b92SHong Zhang +  Arhs - The right-hand side matrix
612cda39b92SHong Zhang .  Alhs - The left-hand side matrix
613cda39b92SHong Zhang -  ctx - User-defined context for matrix evaluation routine
614cda39b92SHong Zhang 
615cda39b92SHong Zhang    Notes: You can pass in PETSC_NULL for any return argument you do not need.
616cda39b92SHong Zhang 
617cda39b92SHong Zhang    Level: intermediate
618cda39b92SHong Zhang 
619cda39b92SHong Zhang .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
620cda39b92SHong Zhang 
621cda39b92SHong Zhang .keywords: TS, timestep, get, matrix
622cda39b92SHong Zhang 
623cda39b92SHong Zhang @*/
624cda39b92SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx)
625cda39b92SHong Zhang {
626cda39b92SHong Zhang   PetscFunctionBegin;
627*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
628cda39b92SHong Zhang   if (Arhs) *Arhs = ts->Arhs;
629cda39b92SHong Zhang   if (Alhs) *Alhs = ts->Alhs;
630cda39b92SHong Zhang   if (ctx)  *ctx = ts->jacP;
631cda39b92SHong Zhang   PetscFunctionReturn(0);
632cda39b92SHong Zhang }
633cda39b92SHong Zhang 
634cda39b92SHong Zhang #undef __FUNCT__
6354a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
636d763cef2SBarry Smith /*@C
637d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
638d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
63976f2fa84SHong Zhang    Use TSSetMatrices() for linear problems.
640d763cef2SBarry Smith 
641d763cef2SBarry Smith    Collective on TS
642d763cef2SBarry Smith 
643d763cef2SBarry Smith    Input Parameters:
644d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
645d763cef2SBarry Smith .  A   - Jacobian matrix
646d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
647d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
648d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
649d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
650d763cef2SBarry Smith 
651d763cef2SBarry Smith    Calling sequence of func:
65287828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
653d763cef2SBarry Smith 
654d763cef2SBarry Smith +  t - current timestep
655d763cef2SBarry Smith .  u - input vector
656d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
657d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
658d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
65994b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
660d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
661d763cef2SBarry Smith 
662d763cef2SBarry Smith    Notes:
66394b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
664d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
665d763cef2SBarry Smith 
666d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
667d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
66856335db2SHong Zhang    completely new matrix structure (not just different matrix elements)
669d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
670d763cef2SBarry Smith    throughout the global iterations.
671d763cef2SBarry Smith 
672d763cef2SBarry Smith    Level: beginner
673d763cef2SBarry Smith 
674d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
675d763cef2SBarry Smith 
676d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
67776f2fa84SHong Zhang           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices()
678d763cef2SBarry Smith 
679d763cef2SBarry Smith @*/
68063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
681d763cef2SBarry Smith {
682d763cef2SBarry Smith   PetscFunctionBegin;
683*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
684*0700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,2);
685*0700a824SBarry Smith   PetscValidHeaderSpecific(B,MAT_CLASSID,3);
686c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
687c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
688d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
68976f2fa84SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()");
690d763cef2SBarry Smith   }
691d763cef2SBarry Smith 
692000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
693d763cef2SBarry Smith   ts->jacP             = ctx;
6948beb423aSHong Zhang   ts->Arhs             = A;
695d763cef2SBarry Smith   ts->B                = B;
696d763cef2SBarry Smith   PetscFunctionReturn(0);
697d763cef2SBarry Smith }
698d763cef2SBarry Smith 
699316643e7SJed Brown 
700316643e7SJed Brown #undef __FUNCT__
701316643e7SJed Brown #define __FUNCT__ "TSSetIFunction"
702316643e7SJed Brown /*@C
703316643e7SJed Brown    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
704316643e7SJed Brown 
705316643e7SJed Brown    Collective on TS
706316643e7SJed Brown 
707316643e7SJed Brown    Input Parameters:
708316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
709316643e7SJed Brown .  f   - the function evaluation routine
710316643e7SJed Brown -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
711316643e7SJed Brown 
712316643e7SJed Brown    Calling sequence of f:
713316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
714316643e7SJed Brown 
715316643e7SJed Brown +  t   - time at step/stage being solved
716316643e7SJed Brown .  u   - state vector
717316643e7SJed Brown .  u_t - time derivative of state vector
718316643e7SJed Brown .  F   - function vector
719316643e7SJed Brown -  ctx - [optional] user-defined context for matrix evaluation routine
720316643e7SJed Brown 
721316643e7SJed Brown    Important:
722316643e7SJed Brown    The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices().  This routine must be used when not solving an ODE.
723316643e7SJed Brown 
724316643e7SJed Brown    Level: beginner
725316643e7SJed Brown 
726316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian
727316643e7SJed Brown 
728316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian()
729316643e7SJed Brown @*/
730316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIFunction(TS ts,TSIFunction f,void *ctx)
731316643e7SJed Brown {
732316643e7SJed Brown 
733316643e7SJed Brown   PetscFunctionBegin;
734*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
735316643e7SJed Brown   ts->ops->ifunction = f;
736316643e7SJed Brown   ts->funP           = ctx;
737316643e7SJed Brown   PetscFunctionReturn(0);
738316643e7SJed Brown }
739316643e7SJed Brown 
740316643e7SJed Brown #undef __FUNCT__
741316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian"
742316643e7SJed Brown /*@C
743316643e7SJed Brown    TSSetIJacobian - Set the function to compute the Jacobian of
744316643e7SJed Brown    G(U) = F(t,U,U0+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
745316643e7SJed Brown 
746316643e7SJed Brown    Collective on TS
747316643e7SJed Brown 
748316643e7SJed Brown    Input Parameters:
749316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
750316643e7SJed Brown .  A   - Jacobian matrix
751316643e7SJed Brown .  B   - preconditioning matrix for A (may be same as A)
752316643e7SJed Brown .  f   - the Jacobian evaluation routine
753316643e7SJed Brown -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
754316643e7SJed Brown 
755316643e7SJed Brown    Calling sequence of f:
756316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
757316643e7SJed Brown 
758316643e7SJed Brown +  t    - time at step/stage being solved
759316643e7SJed Brown .  u    - state vector
760316643e7SJed Brown .  u_t  - time derivative of state vector
761316643e7SJed Brown .  a    - shift
762316643e7SJed Brown .  A    - Jacobian of G(U) = F(t,U,U0+a*U), equivalent to dF/dU + a*dF/dU_t
763316643e7SJed Brown .  B    - preconditioning matrix for A, may be same as A
764316643e7SJed Brown .  flag - flag indicating information about the preconditioner matrix
765316643e7SJed Brown           structure (same as flag in KSPSetOperators())
766316643e7SJed Brown -  ctx  - [optional] user-defined context for matrix evaluation routine
767316643e7SJed Brown 
768316643e7SJed Brown    Notes:
769316643e7SJed Brown    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
770316643e7SJed Brown 
771316643e7SJed Brown    Level: beginner
772316643e7SJed Brown 
773316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian
774316643e7SJed Brown 
775316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian()
776316643e7SJed Brown 
777316643e7SJed Brown @*/
778316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
779316643e7SJed Brown {
780316643e7SJed Brown   PetscErrorCode ierr;
781316643e7SJed Brown 
782316643e7SJed Brown   PetscFunctionBegin;
783*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
784*0700a824SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
785*0700a824SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
786316643e7SJed Brown   if (A) PetscCheckSameComm(ts,1,A,2);
787316643e7SJed Brown   if (B) PetscCheckSameComm(ts,1,B,3);
788316643e7SJed Brown   if (f)   ts->ops->ijacobian = f;
789316643e7SJed Brown   if (ctx) ts->jacP             = ctx;
790316643e7SJed Brown   if (A) {
791316643e7SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
792316643e7SJed Brown     if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);}
793316643e7SJed Brown     ts->A = A;
794316643e7SJed Brown   }
795dc3f620dSJed Brown #if 0
796dc3f620dSJed Brown   /* The sane and consistent alternative */
797316643e7SJed Brown   if (B) {
798316643e7SJed Brown     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
799316643e7SJed Brown     if (ts->B) {ierr = MatDestroy(ts->B);CHKERRQ(ierr);}
800316643e7SJed Brown     ts->B = B;
801316643e7SJed Brown   }
802dc3f620dSJed Brown #else
803dc3f620dSJed Brown   /* Don't reference B because TSDestroy() doesn't destroy it.  These ownership semantics are awkward and inconsistent. */
804dc3f620dSJed Brown   if (B) ts->B = B;
805dc3f620dSJed Brown #endif
806316643e7SJed Brown   PetscFunctionReturn(0);
807316643e7SJed Brown }
808316643e7SJed Brown 
8094a2ae208SSatish Balay #undef __FUNCT__
8104a2ae208SSatish Balay #define __FUNCT__ "TSView"
8117e2c5f70SBarry Smith /*@C
812d763cef2SBarry Smith     TSView - Prints the TS data structure.
813d763cef2SBarry Smith 
8144c49b128SBarry Smith     Collective on TS
815d763cef2SBarry Smith 
816d763cef2SBarry Smith     Input Parameters:
817d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
818d763cef2SBarry Smith -   viewer - visualization context
819d763cef2SBarry Smith 
820d763cef2SBarry Smith     Options Database Key:
821d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
822d763cef2SBarry Smith 
823d763cef2SBarry Smith     Notes:
824d763cef2SBarry Smith     The available visualization contexts include
825b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
826b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
827d763cef2SBarry Smith          output where only the first processor opens
828d763cef2SBarry Smith          the file.  All other processors send their
829d763cef2SBarry Smith          data to the first processor to print.
830d763cef2SBarry Smith 
831d763cef2SBarry Smith     The user can open an alternative visualization context with
832b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
833d763cef2SBarry Smith 
834d763cef2SBarry Smith     Level: beginner
835d763cef2SBarry Smith 
836d763cef2SBarry Smith .keywords: TS, timestep, view
837d763cef2SBarry Smith 
838b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
839d763cef2SBarry Smith @*/
84063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer)
841d763cef2SBarry Smith {
842dfbe8321SBarry Smith   PetscErrorCode ierr;
843a313700dSBarry Smith   const TSType   type;
84432077d6dSBarry Smith   PetscTruth     iascii,isstring;
845d763cef2SBarry Smith 
846d763cef2SBarry Smith   PetscFunctionBegin;
847*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8483050cee2SBarry Smith   if (!viewer) {
8497adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
8503050cee2SBarry Smith   }
851*0700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
852c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
853fd16b177SBarry Smith 
85432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
855b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
85632077d6dSBarry Smith   if (iascii) {
857b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
858a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
859454a90a3SBarry Smith     if (type) {
860b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
861184914b5SBarry Smith     } else {
862b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
863184914b5SBarry Smith     }
864000e7ae3SMatthew Knepley     if (ts->ops->view) {
865b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
866000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
867b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
868d763cef2SBarry Smith     }
86977431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
870a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
871d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
87277431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
873d763cef2SBarry Smith     }
87477431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
8750f5bd95cSBarry Smith   } else if (isstring) {
876a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
877b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
878d763cef2SBarry Smith   }
879b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
88094b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
881d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
882b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
883d763cef2SBarry Smith   PetscFunctionReturn(0);
884d763cef2SBarry Smith }
885d763cef2SBarry Smith 
886d763cef2SBarry Smith 
8874a2ae208SSatish Balay #undef __FUNCT__
8884a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
889d763cef2SBarry Smith /*@C
890d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
891d763cef2SBarry Smith    the timesteppers.
892d763cef2SBarry Smith 
893d763cef2SBarry Smith    Collective on TS
894d763cef2SBarry Smith 
895d763cef2SBarry Smith    Input Parameters:
896d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
897d763cef2SBarry Smith -  usrP - optional user context
898d763cef2SBarry Smith 
899d763cef2SBarry Smith    Level: intermediate
900d763cef2SBarry Smith 
901d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
902d763cef2SBarry Smith 
903d763cef2SBarry Smith .seealso: TSGetApplicationContext()
904d763cef2SBarry Smith @*/
90563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP)
906d763cef2SBarry Smith {
907d763cef2SBarry Smith   PetscFunctionBegin;
908*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
909d763cef2SBarry Smith   ts->user = usrP;
910d763cef2SBarry Smith   PetscFunctionReturn(0);
911d763cef2SBarry Smith }
912d763cef2SBarry Smith 
9134a2ae208SSatish Balay #undef __FUNCT__
9144a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
915d763cef2SBarry Smith /*@C
916d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
917d763cef2SBarry Smith     timestepper.
918d763cef2SBarry Smith 
919d763cef2SBarry Smith     Not Collective
920d763cef2SBarry Smith 
921d763cef2SBarry Smith     Input Parameter:
922d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
923d763cef2SBarry Smith 
924d763cef2SBarry Smith     Output Parameter:
925d763cef2SBarry Smith .   usrP - user context
926d763cef2SBarry Smith 
927d763cef2SBarry Smith     Level: intermediate
928d763cef2SBarry Smith 
929d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
930d763cef2SBarry Smith 
931d763cef2SBarry Smith .seealso: TSSetApplicationContext()
932d763cef2SBarry Smith @*/
93363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP)
934d763cef2SBarry Smith {
935d763cef2SBarry Smith   PetscFunctionBegin;
936*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
937d763cef2SBarry Smith   *usrP = ts->user;
938d763cef2SBarry Smith   PetscFunctionReturn(0);
939d763cef2SBarry Smith }
940d763cef2SBarry Smith 
9414a2ae208SSatish Balay #undef __FUNCT__
9424a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
943d763cef2SBarry Smith /*@
944d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
945d763cef2SBarry Smith 
946d763cef2SBarry Smith    Not Collective
947d763cef2SBarry Smith 
948d763cef2SBarry Smith    Input Parameter:
949d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
950d763cef2SBarry Smith 
951d763cef2SBarry Smith    Output Parameter:
952d763cef2SBarry Smith .  iter - number steps so far
953d763cef2SBarry Smith 
954d763cef2SBarry Smith    Level: intermediate
955d763cef2SBarry Smith 
956d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
957d763cef2SBarry Smith @*/
95863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter)
959d763cef2SBarry Smith {
960d763cef2SBarry Smith   PetscFunctionBegin;
961*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
9624482741eSBarry Smith   PetscValidIntPointer(iter,2);
963d763cef2SBarry Smith   *iter = ts->steps;
964d763cef2SBarry Smith   PetscFunctionReturn(0);
965d763cef2SBarry Smith }
966d763cef2SBarry Smith 
9674a2ae208SSatish Balay #undef __FUNCT__
9684a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
969d763cef2SBarry Smith /*@
970d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
971d763cef2SBarry Smith    as well as the initial time.
972d763cef2SBarry Smith 
973d763cef2SBarry Smith    Collective on TS
974d763cef2SBarry Smith 
975d763cef2SBarry Smith    Input Parameters:
976d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
977d763cef2SBarry Smith .  initial_time - the initial time
978d763cef2SBarry Smith -  time_step - the size of the timestep
979d763cef2SBarry Smith 
980d763cef2SBarry Smith    Level: intermediate
981d763cef2SBarry Smith 
982d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
983d763cef2SBarry Smith 
984d763cef2SBarry Smith .keywords: TS, set, initial, timestep
985d763cef2SBarry Smith @*/
98663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
987d763cef2SBarry Smith {
988d763cef2SBarry Smith   PetscFunctionBegin;
989*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
990d763cef2SBarry Smith   ts->time_step         = time_step;
991d763cef2SBarry Smith   ts->initial_time_step = time_step;
992d763cef2SBarry Smith   ts->ptime             = initial_time;
993d763cef2SBarry Smith   PetscFunctionReturn(0);
994d763cef2SBarry Smith }
995d763cef2SBarry Smith 
9964a2ae208SSatish Balay #undef __FUNCT__
9974a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
998d763cef2SBarry Smith /*@
999d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
1000d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
1001d763cef2SBarry Smith 
1002d763cef2SBarry Smith    Collective on TS
1003d763cef2SBarry Smith 
1004d763cef2SBarry Smith    Input Parameters:
1005d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1006d763cef2SBarry Smith -  time_step - the size of the timestep
1007d763cef2SBarry Smith 
1008d763cef2SBarry Smith    Level: intermediate
1009d763cef2SBarry Smith 
1010d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1011d763cef2SBarry Smith 
1012d763cef2SBarry Smith .keywords: TS, set, timestep
1013d763cef2SBarry Smith @*/
101463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step)
1015d763cef2SBarry Smith {
1016d763cef2SBarry Smith   PetscFunctionBegin;
1017*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1018d763cef2SBarry Smith   ts->time_step = time_step;
1019d763cef2SBarry Smith   PetscFunctionReturn(0);
1020d763cef2SBarry Smith }
1021d763cef2SBarry Smith 
10224a2ae208SSatish Balay #undef __FUNCT__
10234a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
1024d763cef2SBarry Smith /*@
1025d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
1026d763cef2SBarry Smith 
1027d763cef2SBarry Smith    Not Collective
1028d763cef2SBarry Smith 
1029d763cef2SBarry Smith    Input Parameter:
1030d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1031d763cef2SBarry Smith 
1032d763cef2SBarry Smith    Output Parameter:
1033d763cef2SBarry Smith .  dt - the current timestep size
1034d763cef2SBarry Smith 
1035d763cef2SBarry Smith    Level: intermediate
1036d763cef2SBarry Smith 
1037d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1038d763cef2SBarry Smith 
1039d763cef2SBarry Smith .keywords: TS, get, timestep
1040d763cef2SBarry Smith @*/
104163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt)
1042d763cef2SBarry Smith {
1043d763cef2SBarry Smith   PetscFunctionBegin;
1044*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10454482741eSBarry Smith   PetscValidDoublePointer(dt,2);
1046d763cef2SBarry Smith   *dt = ts->time_step;
1047d763cef2SBarry Smith   PetscFunctionReturn(0);
1048d763cef2SBarry Smith }
1049d763cef2SBarry Smith 
10504a2ae208SSatish Balay #undef __FUNCT__
10514a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
1052d8e5e3e6SSatish Balay /*@
1053d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
1054d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
1055d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
1056d763cef2SBarry Smith    the solution at the next timestep has been calculated.
1057d763cef2SBarry Smith 
1058d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
1059d763cef2SBarry Smith 
1060d763cef2SBarry Smith    Input Parameter:
1061d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1062d763cef2SBarry Smith 
1063d763cef2SBarry Smith    Output Parameter:
1064d763cef2SBarry Smith .  v - the vector containing the solution
1065d763cef2SBarry Smith 
1066d763cef2SBarry Smith    Level: intermediate
1067d763cef2SBarry Smith 
1068d763cef2SBarry Smith .seealso: TSGetTimeStep()
1069d763cef2SBarry Smith 
1070d763cef2SBarry Smith .keywords: TS, timestep, get, solution
1071d763cef2SBarry Smith @*/
107263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v)
1073d763cef2SBarry Smith {
1074d763cef2SBarry Smith   PetscFunctionBegin;
1075*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10764482741eSBarry Smith   PetscValidPointer(v,2);
1077d763cef2SBarry Smith   *v = ts->vec_sol_always;
1078d763cef2SBarry Smith   PetscFunctionReturn(0);
1079d763cef2SBarry Smith }
1080d763cef2SBarry Smith 
1081bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
10824a2ae208SSatish Balay #undef __FUNCT__
1083bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
1084d8e5e3e6SSatish Balay /*@
1085bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
1086d763cef2SBarry Smith 
1087bdad233fSMatthew Knepley   Not collective
1088d763cef2SBarry Smith 
1089bdad233fSMatthew Knepley   Input Parameters:
1090bdad233fSMatthew Knepley + ts   - The TS
1091bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1092d763cef2SBarry Smith .vb
1093d763cef2SBarry Smith          U_t = A U
1094d763cef2SBarry Smith          U_t = A(t) U
1095d763cef2SBarry Smith          U_t = F(t,U)
1096d763cef2SBarry Smith .ve
1097d763cef2SBarry Smith 
1098d763cef2SBarry Smith    Level: beginner
1099d763cef2SBarry Smith 
1100bdad233fSMatthew Knepley .keywords: TS, problem type
1101bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1102d763cef2SBarry Smith @*/
110363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type)
1104a7cc72afSBarry Smith {
1105d763cef2SBarry Smith   PetscFunctionBegin;
1106*0700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1107bdad233fSMatthew Knepley   ts->problem_type = type;
1108d763cef2SBarry Smith   PetscFunctionReturn(0);
1109d763cef2SBarry Smith }
1110d763cef2SBarry Smith 
1111bdad233fSMatthew Knepley #undef __FUNCT__
1112bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
1113bdad233fSMatthew Knepley /*@C
1114bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
1115bdad233fSMatthew Knepley 
1116bdad233fSMatthew Knepley   Not collective
1117bdad233fSMatthew Knepley 
1118bdad233fSMatthew Knepley   Input Parameter:
1119bdad233fSMatthew Knepley . ts   - The TS
1120bdad233fSMatthew Knepley 
1121bdad233fSMatthew Knepley   Output Parameter:
1122bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1123bdad233fSMatthew Knepley .vb
1124bdad233fSMatthew Knepley          U_t = A U
1125bdad233fSMatthew Knepley          U_t = A(t) U
1126bdad233fSMatthew Knepley          U_t = F(t,U)
1127bdad233fSMatthew Knepley .ve
1128bdad233fSMatthew Knepley 
1129bdad233fSMatthew Knepley    Level: beginner
1130bdad233fSMatthew Knepley 
1131bdad233fSMatthew Knepley .keywords: TS, problem type
1132bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1133bdad233fSMatthew Knepley @*/
113463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type)
1135a7cc72afSBarry Smith {
1136bdad233fSMatthew Knepley   PetscFunctionBegin;
1137*0700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
11384482741eSBarry Smith   PetscValidIntPointer(type,2);
1139bdad233fSMatthew Knepley   *type = ts->problem_type;
1140bdad233fSMatthew Knepley   PetscFunctionReturn(0);
1141bdad233fSMatthew Knepley }
1142d763cef2SBarry Smith 
11434a2ae208SSatish Balay #undef __FUNCT__
11444a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
1145d763cef2SBarry Smith /*@
1146d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
1147d763cef2SBarry Smith    of a timestepper.
1148d763cef2SBarry Smith 
1149d763cef2SBarry Smith    Collective on TS
1150d763cef2SBarry Smith 
1151d763cef2SBarry Smith    Input Parameter:
1152d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1153d763cef2SBarry Smith 
1154d763cef2SBarry Smith    Notes:
1155d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
1156d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
1157d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
1158d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
1159d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
1160d763cef2SBarry Smith 
1161d763cef2SBarry Smith    Level: advanced
1162d763cef2SBarry Smith 
1163d763cef2SBarry Smith .keywords: TS, timestep, setup
1164d763cef2SBarry Smith 
1165d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
1166d763cef2SBarry Smith @*/
116763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts)
1168d763cef2SBarry Smith {
1169dfbe8321SBarry Smith   PetscErrorCode ierr;
1170d763cef2SBarry Smith 
1171d763cef2SBarry Smith   PetscFunctionBegin;
1172*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
117329bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
11747adad957SLisandro Dalcin   if (!((PetscObject)ts)->type_name) {
11759596e0b4SJed Brown     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1176d763cef2SBarry Smith   }
1177000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1178d763cef2SBarry Smith   ts->setupcalled = 1;
1179d763cef2SBarry Smith   PetscFunctionReturn(0);
1180d763cef2SBarry Smith }
1181d763cef2SBarry Smith 
11824a2ae208SSatish Balay #undef __FUNCT__
11834a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
1184d8e5e3e6SSatish Balay /*@
1185d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
1186d763cef2SBarry Smith    with TSCreate().
1187d763cef2SBarry Smith 
1188d763cef2SBarry Smith    Collective on TS
1189d763cef2SBarry Smith 
1190d763cef2SBarry Smith    Input Parameter:
1191d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1192d763cef2SBarry Smith 
1193d763cef2SBarry Smith    Level: beginner
1194d763cef2SBarry Smith 
1195d763cef2SBarry Smith .keywords: TS, timestepper, destroy
1196d763cef2SBarry Smith 
1197d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
1198d763cef2SBarry Smith @*/
119963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts)
1200d763cef2SBarry Smith {
12016849ba73SBarry Smith   PetscErrorCode ierr;
1202d763cef2SBarry Smith 
1203d763cef2SBarry Smith   PetscFunctionBegin;
1204*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12057adad957SLisandro Dalcin   if (--((PetscObject)ts)->refct > 0) PetscFunctionReturn(0);
1206d763cef2SBarry Smith 
1207be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
12080f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
12096c699258SBarry Smith 
12106c699258SBarry Smith   if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);}
12118beb423aSHong Zhang   if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr)}
121294b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
1213d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
12141e3347e8SBarry Smith   if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);}
1215a6570f20SBarry Smith   ierr = TSMonitorCancel(ts);CHKERRQ(ierr);
1216a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1217d763cef2SBarry Smith   PetscFunctionReturn(0);
1218d763cef2SBarry Smith }
1219d763cef2SBarry Smith 
12204a2ae208SSatish Balay #undef __FUNCT__
12214a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
1222d8e5e3e6SSatish Balay /*@
1223d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1224d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
1225d763cef2SBarry Smith 
1226d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
1227d763cef2SBarry Smith 
1228d763cef2SBarry Smith    Input Parameter:
1229d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1230d763cef2SBarry Smith 
1231d763cef2SBarry Smith    Output Parameter:
1232d763cef2SBarry Smith .  snes - the nonlinear solver context
1233d763cef2SBarry Smith 
1234d763cef2SBarry Smith    Notes:
1235d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
1236d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
123794b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
1238d763cef2SBarry Smith 
1239d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
1240d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
1241d763cef2SBarry Smith 
1242d763cef2SBarry Smith    Level: beginner
1243d763cef2SBarry Smith 
1244d763cef2SBarry Smith .keywords: timestep, get, SNES
1245d763cef2SBarry Smith @*/
124663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes)
1247d763cef2SBarry Smith {
1248d763cef2SBarry Smith   PetscFunctionBegin;
1249*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12504482741eSBarry Smith   PetscValidPointer(snes,2);
1251447600ffSHong Zhang   if (((PetscObject)ts)->type_name == PETSC_NULL)
1252447600ffSHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first");
125394b7f48cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
1254d763cef2SBarry Smith   *snes = ts->snes;
1255d763cef2SBarry Smith   PetscFunctionReturn(0);
1256d763cef2SBarry Smith }
1257d763cef2SBarry Smith 
12584a2ae208SSatish Balay #undef __FUNCT__
125994b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
1260d8e5e3e6SSatish Balay /*@
126194b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
1262d763cef2SBarry Smith    a TS (timestepper) context.
1263d763cef2SBarry Smith 
126494b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1265d763cef2SBarry Smith 
1266d763cef2SBarry Smith    Input Parameter:
1267d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1268d763cef2SBarry Smith 
1269d763cef2SBarry Smith    Output Parameter:
127094b7f48cSBarry Smith .  ksp - the nonlinear solver context
1271d763cef2SBarry Smith 
1272d763cef2SBarry Smith    Notes:
127394b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1274d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1275d763cef2SBarry Smith    KSP and PC contexts as well.
1276d763cef2SBarry Smith 
127794b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
127894b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1279d763cef2SBarry Smith 
1280d763cef2SBarry Smith    Level: beginner
1281d763cef2SBarry Smith 
128294b7f48cSBarry Smith .keywords: timestep, get, KSP
1283d763cef2SBarry Smith @*/
128463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp)
1285d763cef2SBarry Smith {
1286d763cef2SBarry Smith   PetscFunctionBegin;
1287*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12884482741eSBarry Smith   PetscValidPointer(ksp,2);
1289988402f6SHong Zhang   if (((PetscObject)ts)->type_name == PETSC_NULL)
1290988402f6SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
129129bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
129294b7f48cSBarry Smith   *ksp = ts->ksp;
1293d763cef2SBarry Smith   PetscFunctionReturn(0);
1294d763cef2SBarry Smith }
1295d763cef2SBarry Smith 
1296d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1297d763cef2SBarry Smith 
12984a2ae208SSatish Balay #undef __FUNCT__
1299adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1300adb62b0dSMatthew Knepley /*@
1301adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1302adb62b0dSMatthew Knepley    maximum time for iteration.
1303adb62b0dSMatthew Knepley 
1304adb62b0dSMatthew Knepley    Collective on TS
1305adb62b0dSMatthew Knepley 
1306adb62b0dSMatthew Knepley    Input Parameters:
1307adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1308adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1309adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1310adb62b0dSMatthew Knepley 
1311adb62b0dSMatthew Knepley    Level: intermediate
1312adb62b0dSMatthew Knepley 
1313adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1314adb62b0dSMatthew Knepley @*/
131563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1316adb62b0dSMatthew Knepley {
1317adb62b0dSMatthew Knepley   PetscFunctionBegin;
1318*0700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1319abc0a331SBarry Smith   if (maxsteps) {
13204482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1321adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1322adb62b0dSMatthew Knepley   }
1323abc0a331SBarry Smith   if (maxtime ) {
13244482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1325adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1326adb62b0dSMatthew Knepley   }
1327adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1328adb62b0dSMatthew Knepley }
1329adb62b0dSMatthew Knepley 
1330adb62b0dSMatthew Knepley #undef __FUNCT__
13314a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1332d763cef2SBarry Smith /*@
1333d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1334d763cef2SBarry Smith    maximum time for iteration.
1335d763cef2SBarry Smith 
1336d763cef2SBarry Smith    Collective on TS
1337d763cef2SBarry Smith 
1338d763cef2SBarry Smith    Input Parameters:
1339d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1340d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1341d763cef2SBarry Smith -  maxtime - final time to iterate to
1342d763cef2SBarry Smith 
1343d763cef2SBarry Smith    Options Database Keys:
1344d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1345d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1346d763cef2SBarry Smith 
1347d763cef2SBarry Smith    Notes:
1348d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1349d763cef2SBarry Smith 
1350d763cef2SBarry Smith    Level: intermediate
1351d763cef2SBarry Smith 
1352d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1353d763cef2SBarry Smith @*/
135463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1355d763cef2SBarry Smith {
1356d763cef2SBarry Smith   PetscFunctionBegin;
1357*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1358d763cef2SBarry Smith   ts->max_steps = maxsteps;
1359d763cef2SBarry Smith   ts->max_time  = maxtime;
1360d763cef2SBarry Smith   PetscFunctionReturn(0);
1361d763cef2SBarry Smith }
1362d763cef2SBarry Smith 
13634a2ae208SSatish Balay #undef __FUNCT__
13644a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1365d763cef2SBarry Smith /*@
1366d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1367d763cef2SBarry Smith    for use by the TS routines.
1368d763cef2SBarry Smith 
1369d763cef2SBarry Smith    Collective on TS and Vec
1370d763cef2SBarry Smith 
1371d763cef2SBarry Smith    Input Parameters:
1372d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1373d763cef2SBarry Smith -  x - the solution vector
1374d763cef2SBarry Smith 
1375d763cef2SBarry Smith    Level: beginner
1376d763cef2SBarry Smith 
1377d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1378d763cef2SBarry Smith @*/
137963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x)
1380d763cef2SBarry Smith {
1381d763cef2SBarry Smith   PetscFunctionBegin;
1382*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1383*0700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1384d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1385d763cef2SBarry Smith   PetscFunctionReturn(0);
1386d763cef2SBarry Smith }
1387d763cef2SBarry Smith 
1388e74ef692SMatthew Knepley #undef __FUNCT__
1389e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1390ac226902SBarry Smith /*@C
1391000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
13923f2090d5SJed Brown   called once at the beginning of each time step.
1393000e7ae3SMatthew Knepley 
1394000e7ae3SMatthew Knepley   Collective on TS
1395000e7ae3SMatthew Knepley 
1396000e7ae3SMatthew Knepley   Input Parameters:
1397000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1398000e7ae3SMatthew Knepley - func - The function
1399000e7ae3SMatthew Knepley 
1400000e7ae3SMatthew Knepley   Calling sequence of func:
1401000e7ae3SMatthew Knepley . func (TS ts);
1402000e7ae3SMatthew Knepley 
1403000e7ae3SMatthew Knepley   Level: intermediate
1404000e7ae3SMatthew Knepley 
1405000e7ae3SMatthew Knepley .keywords: TS, timestep
1406000e7ae3SMatthew Knepley @*/
140763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1408000e7ae3SMatthew Knepley {
1409000e7ae3SMatthew Knepley   PetscFunctionBegin;
1410*0700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1411000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1412000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1413000e7ae3SMatthew Knepley }
1414000e7ae3SMatthew Knepley 
1415e74ef692SMatthew Knepley #undef __FUNCT__
14163f2090d5SJed Brown #define __FUNCT__ "TSPreStep"
14173f2090d5SJed Brown /*@C
14183f2090d5SJed Brown   TSPreStep - Runs the user-defined pre-step function.
14193f2090d5SJed Brown 
14203f2090d5SJed Brown   Collective on TS
14213f2090d5SJed Brown 
14223f2090d5SJed Brown   Input Parameters:
14233f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
14243f2090d5SJed Brown 
14253f2090d5SJed Brown   Notes:
14263f2090d5SJed Brown   TSPreStep() is typically used within time stepping implementations,
14273f2090d5SJed Brown   so most users would not generally call this routine themselves.
14283f2090d5SJed Brown 
14293f2090d5SJed Brown   Level: developer
14303f2090d5SJed Brown 
14313f2090d5SJed Brown .keywords: TS, timestep
14323f2090d5SJed Brown @*/
14333f2090d5SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSPreStep(TS ts)
14343f2090d5SJed Brown {
14353f2090d5SJed Brown   PetscErrorCode ierr;
14363f2090d5SJed Brown 
14373f2090d5SJed Brown   PetscFunctionBegin;
1438*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
143972ac3e02SJed Brown   if (ts->ops->prestep) {
14403f2090d5SJed Brown     PetscStackPush("TS PreStep function");
14413f2090d5SJed Brown     CHKMEMQ;
14423f2090d5SJed Brown     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
14433f2090d5SJed Brown     CHKMEMQ;
14443f2090d5SJed Brown     PetscStackPop;
1445312ce896SJed Brown   }
14463f2090d5SJed Brown   PetscFunctionReturn(0);
14473f2090d5SJed Brown }
14483f2090d5SJed Brown 
14493f2090d5SJed Brown #undef __FUNCT__
1450e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1451000e7ae3SMatthew Knepley /*@
1452000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1453000e7ae3SMatthew Knepley 
1454000e7ae3SMatthew Knepley   Collective on TS
1455000e7ae3SMatthew Knepley 
1456000e7ae3SMatthew Knepley   Input Parameters:
1457000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1458000e7ae3SMatthew Knepley 
1459000e7ae3SMatthew Knepley   Level: developer
1460000e7ae3SMatthew Knepley 
1461000e7ae3SMatthew Knepley .keywords: TS, timestep
1462000e7ae3SMatthew Knepley @*/
146363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts)
1464000e7ae3SMatthew Knepley {
1465000e7ae3SMatthew Knepley   PetscFunctionBegin;
1466000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1467000e7ae3SMatthew Knepley }
1468000e7ae3SMatthew Knepley 
1469e74ef692SMatthew Knepley #undef __FUNCT__
1470e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1471ac226902SBarry Smith /*@C
1472000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
14733f2090d5SJed Brown   called once at the end of each time step.
1474000e7ae3SMatthew Knepley 
1475000e7ae3SMatthew Knepley   Collective on TS
1476000e7ae3SMatthew Knepley 
1477000e7ae3SMatthew Knepley   Input Parameters:
1478000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1479000e7ae3SMatthew Knepley - func - The function
1480000e7ae3SMatthew Knepley 
1481000e7ae3SMatthew Knepley   Calling sequence of func:
1482000e7ae3SMatthew Knepley . func (TS ts);
1483000e7ae3SMatthew Knepley 
1484000e7ae3SMatthew Knepley   Level: intermediate
1485000e7ae3SMatthew Knepley 
1486000e7ae3SMatthew Knepley .keywords: TS, timestep
1487000e7ae3SMatthew Knepley @*/
148863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1489000e7ae3SMatthew Knepley {
1490000e7ae3SMatthew Knepley   PetscFunctionBegin;
1491*0700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1492000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1493000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1494000e7ae3SMatthew Knepley }
1495000e7ae3SMatthew Knepley 
1496e74ef692SMatthew Knepley #undef __FUNCT__
14973f2090d5SJed Brown #define __FUNCT__ "TSPostStep"
14983f2090d5SJed Brown /*@C
14993f2090d5SJed Brown   TSPostStep - Runs the user-defined post-step function.
15003f2090d5SJed Brown 
15013f2090d5SJed Brown   Collective on TS
15023f2090d5SJed Brown 
15033f2090d5SJed Brown   Input Parameters:
15043f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
15053f2090d5SJed Brown 
15063f2090d5SJed Brown   Notes:
15073f2090d5SJed Brown   TSPostStep() is typically used within time stepping implementations,
15083f2090d5SJed Brown   so most users would not generally call this routine themselves.
15093f2090d5SJed Brown 
15103f2090d5SJed Brown   Level: developer
15113f2090d5SJed Brown 
15123f2090d5SJed Brown .keywords: TS, timestep
15133f2090d5SJed Brown @*/
15143f2090d5SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSPostStep(TS ts)
15153f2090d5SJed Brown {
15163f2090d5SJed Brown   PetscErrorCode ierr;
15173f2090d5SJed Brown 
15183f2090d5SJed Brown   PetscFunctionBegin;
1519*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
152072ac3e02SJed Brown   if (ts->ops->poststep) {
15213f2090d5SJed Brown     PetscStackPush("TS PostStep function");
15223f2090d5SJed Brown     CHKMEMQ;
15233f2090d5SJed Brown     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
15243f2090d5SJed Brown     CHKMEMQ;
15253f2090d5SJed Brown     PetscStackPop;
152672ac3e02SJed Brown   }
15273f2090d5SJed Brown   PetscFunctionReturn(0);
15283f2090d5SJed Brown }
15293f2090d5SJed Brown 
15303f2090d5SJed Brown #undef __FUNCT__
1531e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1532000e7ae3SMatthew Knepley /*@
1533000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1534000e7ae3SMatthew Knepley 
1535000e7ae3SMatthew Knepley   Collective on TS
1536000e7ae3SMatthew Knepley 
1537000e7ae3SMatthew Knepley   Input Parameters:
1538000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1539000e7ae3SMatthew Knepley 
1540000e7ae3SMatthew Knepley   Level: developer
1541000e7ae3SMatthew Knepley 
1542000e7ae3SMatthew Knepley .keywords: TS, timestep
1543000e7ae3SMatthew Knepley @*/
154463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts)
1545000e7ae3SMatthew Knepley {
1546000e7ae3SMatthew Knepley   PetscFunctionBegin;
1547000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1548000e7ae3SMatthew Knepley }
1549000e7ae3SMatthew Knepley 
1550d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1551d763cef2SBarry Smith 
15524a2ae208SSatish Balay #undef __FUNCT__
1553a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet"
1554d763cef2SBarry Smith /*@C
1555a6570f20SBarry Smith    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1556d763cef2SBarry Smith    timestep to display the iteration's  progress.
1557d763cef2SBarry Smith 
1558d763cef2SBarry Smith    Collective on TS
1559d763cef2SBarry Smith 
1560d763cef2SBarry Smith    Input Parameters:
1561d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1562d763cef2SBarry Smith .  func - monitoring routine
1563329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1564b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1565b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1566b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1567d763cef2SBarry Smith 
1568d763cef2SBarry Smith    Calling sequence of func:
1569a7cc72afSBarry Smith $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1570d763cef2SBarry Smith 
1571d763cef2SBarry Smith +    ts - the TS context
1572d763cef2SBarry Smith .    steps - iteration number
15731f06c33eSBarry Smith .    time - current time
1574d763cef2SBarry Smith .    x - current iterate
1575d763cef2SBarry Smith -    mctx - [optional] monitoring context
1576d763cef2SBarry Smith 
1577d763cef2SBarry Smith    Notes:
1578d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1579d763cef2SBarry Smith    already has been loaded.
1580d763cef2SBarry Smith 
1581025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each TS object
1582025f1a04SBarry Smith 
1583d763cef2SBarry Smith    Level: intermediate
1584d763cef2SBarry Smith 
1585d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1586d763cef2SBarry Smith 
1587a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel()
1588d763cef2SBarry Smith @*/
1589a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1590d763cef2SBarry Smith {
1591d763cef2SBarry Smith   PetscFunctionBegin;
1592*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1593d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
159429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1595d763cef2SBarry Smith   }
1596d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1597329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1598d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1599d763cef2SBarry Smith   PetscFunctionReturn(0);
1600d763cef2SBarry Smith }
1601d763cef2SBarry Smith 
16024a2ae208SSatish Balay #undef __FUNCT__
1603a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel"
1604d763cef2SBarry Smith /*@C
1605a6570f20SBarry Smith    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1606d763cef2SBarry Smith 
1607d763cef2SBarry Smith    Collective on TS
1608d763cef2SBarry Smith 
1609d763cef2SBarry Smith    Input Parameters:
1610d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1611d763cef2SBarry Smith 
1612d763cef2SBarry Smith    Notes:
1613d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1614d763cef2SBarry Smith 
1615d763cef2SBarry Smith    Level: intermediate
1616d763cef2SBarry Smith 
1617d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1618d763cef2SBarry Smith 
1619a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet()
1620d763cef2SBarry Smith @*/
1621a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts)
1622d763cef2SBarry Smith {
1623d952e501SBarry Smith   PetscErrorCode ierr;
1624d952e501SBarry Smith   PetscInt       i;
1625d952e501SBarry Smith 
1626d763cef2SBarry Smith   PetscFunctionBegin;
1627*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1628d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1629d952e501SBarry Smith     if (ts->mdestroy[i]) {
1630d952e501SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
1631d952e501SBarry Smith     }
1632d952e501SBarry Smith   }
1633d763cef2SBarry Smith   ts->numbermonitors = 0;
1634d763cef2SBarry Smith   PetscFunctionReturn(0);
1635d763cef2SBarry Smith }
1636d763cef2SBarry Smith 
16374a2ae208SSatish Balay #undef __FUNCT__
1638a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault"
1639d8e5e3e6SSatish Balay /*@
1640a6570f20SBarry Smith    TSMonitorDefault - Sets the Default monitor
16415516499fSSatish Balay 
16425516499fSSatish Balay    Level: intermediate
164341251cbbSSatish Balay 
16445516499fSSatish Balay .keywords: TS, set, monitor
16455516499fSSatish Balay 
164641251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet()
164741251cbbSSatish Balay @*/
1648a6570f20SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1649d763cef2SBarry Smith {
1650dfbe8321SBarry Smith   PetscErrorCode          ierr;
1651a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
1652d132466eSBarry Smith 
1653d763cef2SBarry Smith   PetscFunctionBegin;
1654a34d58ebSBarry Smith   if (!ctx) {
16557adad957SLisandro Dalcin     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
1656a34d58ebSBarry Smith   }
1657a34d58ebSBarry Smith   ierr = PetscViewerASCIIMonitorPrintf(viewer,"timestep %D dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1658a34d58ebSBarry Smith   if (!ctx) {
1659a34d58ebSBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
1660a34d58ebSBarry Smith   }
1661d763cef2SBarry Smith   PetscFunctionReturn(0);
1662d763cef2SBarry Smith }
1663d763cef2SBarry Smith 
16644a2ae208SSatish Balay #undef __FUNCT__
16654a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1666d763cef2SBarry Smith /*@
1667d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1668d763cef2SBarry Smith 
1669d763cef2SBarry Smith    Collective on TS
1670d763cef2SBarry Smith 
1671d763cef2SBarry Smith    Input Parameter:
1672d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1673d763cef2SBarry Smith 
1674d763cef2SBarry Smith    Output Parameters:
1675d763cef2SBarry Smith +  steps - number of iterations until termination
1676142b95e3SSatish Balay -  ptime - time until termination
1677d763cef2SBarry Smith 
1678d763cef2SBarry Smith    Level: beginner
1679d763cef2SBarry Smith 
1680d763cef2SBarry Smith .keywords: TS, timestep, solve
1681d763cef2SBarry Smith 
1682d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1683d763cef2SBarry Smith @*/
168463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1685d763cef2SBarry Smith {
1686dfbe8321SBarry Smith   PetscErrorCode ierr;
1687d763cef2SBarry Smith 
1688d763cef2SBarry Smith   PetscFunctionBegin;
1689*0700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1690d405a339SMatthew Knepley   if (!ts->setupcalled) {
1691d405a339SMatthew Knepley     ierr = TSSetUp(ts);CHKERRQ(ierr);
1692d405a339SMatthew Knepley   }
1693d405a339SMatthew Knepley 
1694d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1695000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1696d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1697d405a339SMatthew Knepley 
16984bb05414SBarry Smith   if (!PetscPreLoadingOn) {
16997adad957SLisandro Dalcin     ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr);
1700d405a339SMatthew Knepley   }
1701d763cef2SBarry Smith   PetscFunctionReturn(0);
1702d763cef2SBarry Smith }
1703d763cef2SBarry Smith 
17044a2ae208SSatish Balay #undef __FUNCT__
17056a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve"
17066a4d4014SLisandro Dalcin /*@
17076a4d4014SLisandro Dalcin    TSSolve - Steps the requested number of timesteps.
17086a4d4014SLisandro Dalcin 
17096a4d4014SLisandro Dalcin    Collective on TS
17106a4d4014SLisandro Dalcin 
17116a4d4014SLisandro Dalcin    Input Parameter:
17126a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
17136a4d4014SLisandro Dalcin -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
17146a4d4014SLisandro Dalcin 
17156a4d4014SLisandro Dalcin    Level: beginner
17166a4d4014SLisandro Dalcin 
17176a4d4014SLisandro Dalcin .keywords: TS, timestep, solve
17186a4d4014SLisandro Dalcin 
17196a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep()
17206a4d4014SLisandro Dalcin @*/
17216a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x)
17226a4d4014SLisandro Dalcin {
17236a4d4014SLisandro Dalcin   PetscInt       steps;
17246a4d4014SLisandro Dalcin   PetscReal      ptime;
17256a4d4014SLisandro Dalcin   PetscErrorCode ierr;
17266a4d4014SLisandro Dalcin   PetscFunctionBegin;
1727*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
17286a4d4014SLisandro Dalcin   /* set solution vector if provided */
17296a4d4014SLisandro Dalcin   if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); }
17306a4d4014SLisandro Dalcin   /* reset time step and iteration counters */
17316a4d4014SLisandro Dalcin   ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0;
17326a4d4014SLisandro Dalcin   /* steps the requested number of timesteps. */
17336a4d4014SLisandro Dalcin   ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr);
17346a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
17356a4d4014SLisandro Dalcin }
17366a4d4014SLisandro Dalcin 
17376a4d4014SLisandro Dalcin #undef __FUNCT__
17384a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1739d763cef2SBarry Smith /*
1740d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1741d763cef2SBarry Smith */
1742a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1743d763cef2SBarry Smith {
17446849ba73SBarry Smith   PetscErrorCode ierr;
1745a7cc72afSBarry Smith   PetscInt i,n = ts->numbermonitors;
1746d763cef2SBarry Smith 
1747d763cef2SBarry Smith   PetscFunctionBegin;
1748d763cef2SBarry Smith   for (i=0; i<n; i++) {
1749142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1750d763cef2SBarry Smith   }
1751d763cef2SBarry Smith   PetscFunctionReturn(0);
1752d763cef2SBarry Smith }
1753d763cef2SBarry Smith 
1754d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1755d763cef2SBarry Smith 
17564a2ae208SSatish Balay #undef __FUNCT__
1757a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate"
1758d763cef2SBarry Smith /*@C
1759a6570f20SBarry Smith    TSMonitorLGCreate - Creates a line graph context for use with
1760d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1761d763cef2SBarry Smith 
1762d763cef2SBarry Smith    Collective on TS
1763d763cef2SBarry Smith 
1764d763cef2SBarry Smith    Input Parameters:
1765d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1766d763cef2SBarry Smith .  label - the title to put in the title bar
17677c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1768d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1769d763cef2SBarry Smith 
1770d763cef2SBarry Smith    Output Parameter:
1771d763cef2SBarry Smith .  draw - the drawing context
1772d763cef2SBarry Smith 
1773d763cef2SBarry Smith    Options Database Key:
1774a6570f20SBarry Smith .  -ts_monitor_draw - automatically sets line graph monitor
1775d763cef2SBarry Smith 
1776d763cef2SBarry Smith    Notes:
1777a6570f20SBarry Smith    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1778d763cef2SBarry Smith 
1779d763cef2SBarry Smith    Level: intermediate
1780d763cef2SBarry Smith 
17817c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1782d763cef2SBarry Smith 
1783a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet()
17847c922b88SBarry Smith 
1785d763cef2SBarry Smith @*/
1786a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1787d763cef2SBarry Smith {
1788b0a32e0cSBarry Smith   PetscDraw      win;
1789dfbe8321SBarry Smith   PetscErrorCode ierr;
1790d763cef2SBarry Smith 
1791d763cef2SBarry Smith   PetscFunctionBegin;
1792b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1793b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1794b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1795b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1796d763cef2SBarry Smith 
179752e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1798d763cef2SBarry Smith   PetscFunctionReturn(0);
1799d763cef2SBarry Smith }
1800d763cef2SBarry Smith 
18014a2ae208SSatish Balay #undef __FUNCT__
1802a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG"
1803a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1804d763cef2SBarry Smith {
1805b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
180687828ca2SBarry Smith   PetscReal      x,y = ptime;
1807dfbe8321SBarry Smith   PetscErrorCode ierr;
1808d763cef2SBarry Smith 
1809d763cef2SBarry Smith   PetscFunctionBegin;
18107c922b88SBarry Smith   if (!monctx) {
18117c922b88SBarry Smith     MPI_Comm    comm;
1812b0a32e0cSBarry Smith     PetscViewer viewer;
18137c922b88SBarry Smith 
18147c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1815b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1816b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
18177c922b88SBarry Smith   }
18187c922b88SBarry Smith 
1819b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
182087828ca2SBarry Smith   x = (PetscReal)n;
1821b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1822d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1823b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1824d763cef2SBarry Smith   }
1825d763cef2SBarry Smith   PetscFunctionReturn(0);
1826d763cef2SBarry Smith }
1827d763cef2SBarry Smith 
18284a2ae208SSatish Balay #undef __FUNCT__
1829a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy"
1830d763cef2SBarry Smith /*@C
1831a6570f20SBarry Smith    TSMonitorLGDestroy - Destroys a line graph context that was created
1832a6570f20SBarry Smith    with TSMonitorLGCreate().
1833d763cef2SBarry Smith 
1834b0a32e0cSBarry Smith    Collective on PetscDrawLG
1835d763cef2SBarry Smith 
1836d763cef2SBarry Smith    Input Parameter:
1837d763cef2SBarry Smith .  draw - the drawing context
1838d763cef2SBarry Smith 
1839d763cef2SBarry Smith    Level: intermediate
1840d763cef2SBarry Smith 
1841d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1842d763cef2SBarry Smith 
1843a6570f20SBarry Smith .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1844d763cef2SBarry Smith @*/
1845a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg)
1846d763cef2SBarry Smith {
1847b0a32e0cSBarry Smith   PetscDraw      draw;
1848dfbe8321SBarry Smith   PetscErrorCode ierr;
1849d763cef2SBarry Smith 
1850d763cef2SBarry Smith   PetscFunctionBegin;
1851b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1852b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1853b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1854d763cef2SBarry Smith   PetscFunctionReturn(0);
1855d763cef2SBarry Smith }
1856d763cef2SBarry Smith 
18574a2ae208SSatish Balay #undef __FUNCT__
18584a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1859d763cef2SBarry Smith /*@
1860d763cef2SBarry Smith    TSGetTime - Gets the current time.
1861d763cef2SBarry Smith 
1862d763cef2SBarry Smith    Not Collective
1863d763cef2SBarry Smith 
1864d763cef2SBarry Smith    Input Parameter:
1865d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1866d763cef2SBarry Smith 
1867d763cef2SBarry Smith    Output Parameter:
1868d763cef2SBarry Smith .  t  - the current time
1869d763cef2SBarry Smith 
1870d763cef2SBarry Smith    Level: beginner
1871d763cef2SBarry Smith 
1872d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1873d763cef2SBarry Smith 
1874d763cef2SBarry Smith .keywords: TS, get, time
1875d763cef2SBarry Smith @*/
187663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t)
1877d763cef2SBarry Smith {
1878d763cef2SBarry Smith   PetscFunctionBegin;
1879*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
18804482741eSBarry Smith   PetscValidDoublePointer(t,2);
1881d763cef2SBarry Smith   *t = ts->ptime;
1882d763cef2SBarry Smith   PetscFunctionReturn(0);
1883d763cef2SBarry Smith }
1884d763cef2SBarry Smith 
18854a2ae208SSatish Balay #undef __FUNCT__
18866a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime"
18876a4d4014SLisandro Dalcin /*@
18886a4d4014SLisandro Dalcin    TSSetTime - Allows one to reset the time.
18896a4d4014SLisandro Dalcin 
18906a4d4014SLisandro Dalcin    Collective on TS
18916a4d4014SLisandro Dalcin 
18926a4d4014SLisandro Dalcin    Input Parameters:
18936a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
18946a4d4014SLisandro Dalcin -  time - the time
18956a4d4014SLisandro Dalcin 
18966a4d4014SLisandro Dalcin    Level: intermediate
18976a4d4014SLisandro Dalcin 
18986a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration()
18996a4d4014SLisandro Dalcin 
19006a4d4014SLisandro Dalcin .keywords: TS, set, time
19016a4d4014SLisandro Dalcin @*/
19026a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t)
19036a4d4014SLisandro Dalcin {
19046a4d4014SLisandro Dalcin   PetscFunctionBegin;
1905*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
19066a4d4014SLisandro Dalcin   ts->ptime = t;
19076a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
19086a4d4014SLisandro Dalcin }
19096a4d4014SLisandro Dalcin 
19106a4d4014SLisandro Dalcin #undef __FUNCT__
19114a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1912d763cef2SBarry Smith /*@C
1913d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1914d763cef2SBarry Smith    TS options in the database.
1915d763cef2SBarry Smith 
1916d763cef2SBarry Smith    Collective on TS
1917d763cef2SBarry Smith 
1918d763cef2SBarry Smith    Input Parameter:
1919d763cef2SBarry Smith +  ts     - The TS context
1920d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1921d763cef2SBarry Smith 
1922d763cef2SBarry Smith    Notes:
1923d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1924d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1925d763cef2SBarry Smith    hyphen.
1926d763cef2SBarry Smith 
1927d763cef2SBarry Smith    Level: advanced
1928d763cef2SBarry Smith 
1929d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1930d763cef2SBarry Smith 
1931d763cef2SBarry Smith .seealso: TSSetFromOptions()
1932d763cef2SBarry Smith 
1933d763cef2SBarry Smith @*/
193463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[])
1935d763cef2SBarry Smith {
1936dfbe8321SBarry Smith   PetscErrorCode ierr;
1937d763cef2SBarry Smith 
1938d763cef2SBarry Smith   PetscFunctionBegin;
1939*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1940d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1941d763cef2SBarry Smith   switch(ts->problem_type) {
1942d763cef2SBarry Smith     case TS_NONLINEAR:
194359580b9cSBarry Smith       if (ts->snes) {
1944d763cef2SBarry Smith         ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
194559580b9cSBarry Smith       }
1946d763cef2SBarry Smith       break;
1947d763cef2SBarry Smith     case TS_LINEAR:
194859580b9cSBarry Smith       if (ts->ksp) {
194994b7f48cSBarry Smith         ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
195059580b9cSBarry Smith       }
1951d763cef2SBarry Smith       break;
1952d763cef2SBarry Smith   }
1953d763cef2SBarry Smith   PetscFunctionReturn(0);
1954d763cef2SBarry Smith }
1955d763cef2SBarry Smith 
1956d763cef2SBarry Smith 
19574a2ae208SSatish Balay #undef __FUNCT__
19584a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1959d763cef2SBarry Smith /*@C
1960d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1961d763cef2SBarry Smith    TS options in the database.
1962d763cef2SBarry Smith 
1963d763cef2SBarry Smith    Collective on TS
1964d763cef2SBarry Smith 
1965d763cef2SBarry Smith    Input Parameter:
1966d763cef2SBarry Smith +  ts     - The TS context
1967d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1968d763cef2SBarry Smith 
1969d763cef2SBarry Smith    Notes:
1970d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1971d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1972d763cef2SBarry Smith    hyphen.
1973d763cef2SBarry Smith 
1974d763cef2SBarry Smith    Level: advanced
1975d763cef2SBarry Smith 
1976d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1977d763cef2SBarry Smith 
1978d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1979d763cef2SBarry Smith 
1980d763cef2SBarry Smith @*/
198163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[])
1982d763cef2SBarry Smith {
1983dfbe8321SBarry Smith   PetscErrorCode ierr;
1984d763cef2SBarry Smith 
1985d763cef2SBarry Smith   PetscFunctionBegin;
1986*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1987d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1988d763cef2SBarry Smith   switch(ts->problem_type) {
1989d763cef2SBarry Smith     case TS_NONLINEAR:
19901ac94b3bSBarry Smith       if (ts->snes) {
1991d763cef2SBarry Smith         ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
19921ac94b3bSBarry Smith       }
1993d763cef2SBarry Smith       break;
1994d763cef2SBarry Smith     case TS_LINEAR:
19951ac94b3bSBarry Smith       if (ts->ksp) {
199694b7f48cSBarry Smith         ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
19971ac94b3bSBarry Smith       }
1998d763cef2SBarry Smith       break;
1999d763cef2SBarry Smith   }
2000d763cef2SBarry Smith   PetscFunctionReturn(0);
2001d763cef2SBarry Smith }
2002d763cef2SBarry Smith 
20034a2ae208SSatish Balay #undef __FUNCT__
20044a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
2005d763cef2SBarry Smith /*@C
2006d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
2007d763cef2SBarry Smith    TS options in the database.
2008d763cef2SBarry Smith 
2009d763cef2SBarry Smith    Not Collective
2010d763cef2SBarry Smith 
2011d763cef2SBarry Smith    Input Parameter:
2012d763cef2SBarry Smith .  ts - The TS context
2013d763cef2SBarry Smith 
2014d763cef2SBarry Smith    Output Parameter:
2015d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
2016d763cef2SBarry Smith 
2017d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
2018d763cef2SBarry Smith    sufficient length to hold the prefix.
2019d763cef2SBarry Smith 
2020d763cef2SBarry Smith    Level: intermediate
2021d763cef2SBarry Smith 
2022d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
2023d763cef2SBarry Smith 
2024d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
2025d763cef2SBarry Smith @*/
2026e060cb09SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[])
2027d763cef2SBarry Smith {
2028dfbe8321SBarry Smith   PetscErrorCode ierr;
2029d763cef2SBarry Smith 
2030d763cef2SBarry Smith   PetscFunctionBegin;
2031*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
20324482741eSBarry Smith   PetscValidPointer(prefix,2);
2033d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2034d763cef2SBarry Smith   PetscFunctionReturn(0);
2035d763cef2SBarry Smith }
2036d763cef2SBarry Smith 
20374a2ae208SSatish Balay #undef __FUNCT__
20384a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
2039d763cef2SBarry Smith /*@C
2040d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2041d763cef2SBarry Smith 
2042d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
2043d763cef2SBarry Smith 
2044d763cef2SBarry Smith    Input Parameter:
2045d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
2046d763cef2SBarry Smith 
2047d763cef2SBarry Smith    Output Parameters:
2048d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
2049d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
2050d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
2051d763cef2SBarry Smith 
2052d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2053d763cef2SBarry Smith 
2054d763cef2SBarry Smith    Level: intermediate
2055d763cef2SBarry Smith 
205626d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2057d763cef2SBarry Smith 
2058d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
2059d763cef2SBarry Smith @*/
206063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
2061d763cef2SBarry Smith {
2062d763cef2SBarry Smith   PetscFunctionBegin;
206326d46c62SHong Zhang   if (J) *J = ts->Arhs;
206426d46c62SHong Zhang   if (M) *M = ts->B;
206526d46c62SHong Zhang   if (ctx) *ctx = ts->jacP;
2066d763cef2SBarry Smith   PetscFunctionReturn(0);
2067d763cef2SBarry Smith }
2068d763cef2SBarry Smith 
20691713a123SBarry Smith #undef __FUNCT__
20702eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian"
20712eca1d9cSJed Brown /*@C
20722eca1d9cSJed Brown    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
20732eca1d9cSJed Brown 
20742eca1d9cSJed Brown    Not Collective, but parallel objects are returned if TS is parallel
20752eca1d9cSJed Brown 
20762eca1d9cSJed Brown    Input Parameter:
20772eca1d9cSJed Brown .  ts  - The TS context obtained from TSCreate()
20782eca1d9cSJed Brown 
20792eca1d9cSJed Brown    Output Parameters:
20802eca1d9cSJed Brown +  A   - The Jacobian of F(t,U,U_t)
20812eca1d9cSJed Brown .  B   - The preconditioner matrix, often the same as A
20822eca1d9cSJed Brown .  f   - The function to compute the matrices
20832eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine
20842eca1d9cSJed Brown 
20852eca1d9cSJed Brown    Notes: You can pass in PETSC_NULL for any return argument you do not need.
20862eca1d9cSJed Brown 
20872eca1d9cSJed Brown    Level: advanced
20882eca1d9cSJed Brown 
20892eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
20902eca1d9cSJed Brown 
20912eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian
20922eca1d9cSJed Brown @*/
20932eca1d9cSJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
20942eca1d9cSJed Brown {
20952eca1d9cSJed Brown   PetscFunctionBegin;
20962eca1d9cSJed Brown   if (A) *A = ts->A;
20972eca1d9cSJed Brown   if (B) *B = ts->B;
20982eca1d9cSJed Brown   if (f) *f = ts->ops->ijacobian;
20992eca1d9cSJed Brown   if (ctx) *ctx = ts->jacP;
21002eca1d9cSJed Brown   PetscFunctionReturn(0);
21012eca1d9cSJed Brown }
21022eca1d9cSJed Brown 
21032eca1d9cSJed Brown #undef __FUNCT__
2104a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution"
21051713a123SBarry Smith /*@C
2106a6570f20SBarry Smith    TSMonitorSolution - Monitors progress of the TS solvers by calling
21071713a123SBarry Smith    VecView() for the solution at each timestep
21081713a123SBarry Smith 
21091713a123SBarry Smith    Collective on TS
21101713a123SBarry Smith 
21111713a123SBarry Smith    Input Parameters:
21121713a123SBarry Smith +  ts - the TS context
21131713a123SBarry Smith .  step - current time-step
2114142b95e3SSatish Balay .  ptime - current time
21151713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
21161713a123SBarry Smith 
21171713a123SBarry Smith    Level: intermediate
21181713a123SBarry Smith 
21191713a123SBarry Smith .keywords: TS,  vector, monitor, view
21201713a123SBarry Smith 
2121a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
21221713a123SBarry Smith @*/
2123a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
21241713a123SBarry Smith {
2125dfbe8321SBarry Smith   PetscErrorCode ierr;
21261713a123SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
21271713a123SBarry Smith 
21281713a123SBarry Smith   PetscFunctionBegin;
2129a34d58ebSBarry Smith   if (!dummy) {
21307adad957SLisandro Dalcin     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
21311713a123SBarry Smith   }
21321713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
21331713a123SBarry Smith   PetscFunctionReturn(0);
21341713a123SBarry Smith }
21351713a123SBarry Smith 
21361713a123SBarry Smith 
21376c699258SBarry Smith #undef __FUNCT__
21386c699258SBarry Smith #define __FUNCT__ "TSSetDM"
21396c699258SBarry Smith /*@
21406c699258SBarry Smith    TSSetDM - Sets the DM that may be used by some preconditioners
21416c699258SBarry Smith 
21426c699258SBarry Smith    Collective on TS
21436c699258SBarry Smith 
21446c699258SBarry Smith    Input Parameters:
21456c699258SBarry Smith +  ts - the preconditioner context
21466c699258SBarry Smith -  dm - the dm
21476c699258SBarry Smith 
21486c699258SBarry Smith    Level: intermediate
21496c699258SBarry Smith 
21506c699258SBarry Smith 
21516c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
21526c699258SBarry Smith @*/
21536c699258SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSSetDM(TS ts,DM dm)
21546c699258SBarry Smith {
21556c699258SBarry Smith   PetscErrorCode ierr;
21566c699258SBarry Smith 
21576c699258SBarry Smith   PetscFunctionBegin;
2158*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
21596c699258SBarry Smith   if (ts->dm) {ierr = DMDestroy(ts->dm);CHKERRQ(ierr);}
21606c699258SBarry Smith   ts->dm = dm;
21616c699258SBarry Smith   ierr = PetscObjectReference((PetscObject)ts->dm);CHKERRQ(ierr);
21626c699258SBarry Smith   if (ts->snes) {ierr = SNESSetDM(ts->snes,dm);CHKERRQ(ierr);}
21636c699258SBarry Smith   PetscFunctionReturn(0);
21646c699258SBarry Smith }
21656c699258SBarry Smith 
21666c699258SBarry Smith #undef __FUNCT__
21676c699258SBarry Smith #define __FUNCT__ "TSGetDM"
21686c699258SBarry Smith /*@
21696c699258SBarry Smith    TSGetDM - Gets the DM that may be used by some preconditioners
21706c699258SBarry Smith 
21716c699258SBarry Smith    Collective on TS
21726c699258SBarry Smith 
21736c699258SBarry Smith    Input Parameter:
21746c699258SBarry Smith . ts - the preconditioner context
21756c699258SBarry Smith 
21766c699258SBarry Smith    Output Parameter:
21776c699258SBarry Smith .  dm - the dm
21786c699258SBarry Smith 
21796c699258SBarry Smith    Level: intermediate
21806c699258SBarry Smith 
21816c699258SBarry Smith 
21826c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
21836c699258SBarry Smith @*/
21846c699258SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetDM(TS ts,DM *dm)
21856c699258SBarry Smith {
21866c699258SBarry Smith   PetscFunctionBegin;
2187*0700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
21886c699258SBarry Smith   *dm = ts->dm;
21896c699258SBarry Smith   PetscFunctionReturn(0);
21906c699258SBarry Smith }
21911713a123SBarry Smith 
2192