xref: /petsc/src/ts/interface/ts.c (revision 63495f91a95b127edaff8973ef60416f316a5d76)
163dd3a1aSKris Buschelman #define PETSCTS_DLL
263dd3a1aSKris Buschelman 
37c4f633dSBarry Smith #include "private/tsimpl.h"        /*I "petscts.h"  I*/
4d763cef2SBarry Smith 
5d5ba7fb7SMatthew Knepley /* Logging support */
6166c7f25SBarry Smith PetscCookie PETSCTS_DLLEXPORT TS_COOKIE;
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 {
35bdad233fSMatthew Knepley     defaultType = TS_EULER;
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:
590f3b3ca1SHong Zhang +  -ts_type <type> - TS_EULER, TS_BEULER, TS_SUNDIALS, TS_PSEUDO, TS_CRANK_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;
814482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,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
1327c236d22SBarry Smith          before they allow users to set options */
1338beb423aSHong Zhang       ierr = SNESSetJacobian(ts->snes,ts->Arhs,ts->B,0,ts);CHKERRQ(ierr);
134bdad233fSMatthew Knepley       ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
135156fc9a6SMatthew Knepley     }
136bdad233fSMatthew Knepley     break;
137bdad233fSMatthew Knepley   default:
13877431f27SBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid problem type: %d", (int)ts->problem_type);
139bdad233fSMatthew Knepley   }
140bdad233fSMatthew Knepley 
141bdad233fSMatthew Knepley   PetscFunctionReturn(0);
142bdad233fSMatthew Knepley }
143bdad233fSMatthew Knepley 
144bdad233fSMatthew Knepley #undef  __FUNCT__
145bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
146bdad233fSMatthew Knepley /*@
147bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
148bdad233fSMatthew Knepley 
149bdad233fSMatthew Knepley   Collective on TS
150bdad233fSMatthew Knepley 
151bdad233fSMatthew Knepley   Input Parameter:
152bdad233fSMatthew Knepley . ts - The ts
153bdad233fSMatthew Knepley 
154bdad233fSMatthew Knepley   Level: intermediate
155bdad233fSMatthew Knepley 
156bdad233fSMatthew Knepley .keywords: TS, view, options, database
157bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
158bdad233fSMatthew Knepley @*/
15963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSViewFromOptions(TS ts,const char title[])
160bdad233fSMatthew Knepley {
161bdad233fSMatthew Knepley   PetscViewer    viewer;
162bdad233fSMatthew Knepley   PetscDraw      draw;
16390d69ab7SBarry Smith   PetscTruth     opt = PETSC_FALSE;
164e10c49a3SBarry Smith   char           fileName[PETSC_MAX_PATH_LEN];
165dfbe8321SBarry Smith   PetscErrorCode ierr;
166bdad233fSMatthew Knepley 
167bdad233fSMatthew Knepley   PetscFunctionBegin;
1687adad957SLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr);
169eabae89aSBarry Smith   if (opt && !PetscPreLoadingOn) {
1707adad957SLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr);
171bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
172bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
173bdad233fSMatthew Knepley   }
1748e83347fSKai Germaschewski   opt = PETSC_FALSE;
17590d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr);
176a7cc72afSBarry Smith   if (opt) {
1777adad957SLisandro Dalcin     ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
178bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
179a7cc72afSBarry Smith     if (title) {
1801836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
181bdad233fSMatthew Knepley     } else {
182bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr);
1837adad957SLisandro Dalcin       ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr);
184bdad233fSMatthew Knepley     }
185bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
186bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
187bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
188bdad233fSMatthew Knepley     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
189bdad233fSMatthew Knepley   }
190bdad233fSMatthew Knepley   PetscFunctionReturn(0);
191bdad233fSMatthew Knepley }
192bdad233fSMatthew Knepley 
193bdad233fSMatthew Knepley #undef __FUNCT__
1944a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
195a7a1495cSBarry Smith /*@
1968c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
197a7a1495cSBarry Smith       set with TSSetRHSJacobian().
198a7a1495cSBarry Smith 
199a7a1495cSBarry Smith    Collective on TS and Vec
200a7a1495cSBarry Smith 
201a7a1495cSBarry Smith    Input Parameters:
202316643e7SJed Brown +  ts - the TS context
203a7a1495cSBarry Smith .  t - current timestep
204a7a1495cSBarry Smith -  x - input vector
205a7a1495cSBarry Smith 
206a7a1495cSBarry Smith    Output Parameters:
207a7a1495cSBarry Smith +  A - Jacobian matrix
208a7a1495cSBarry Smith .  B - optional preconditioning matrix
209a7a1495cSBarry Smith -  flag - flag indicating matrix structure
210a7a1495cSBarry Smith 
211a7a1495cSBarry Smith    Notes:
212a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
213a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
214a7a1495cSBarry Smith 
21594b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
216a7a1495cSBarry Smith    flag parameter.
217a7a1495cSBarry Smith 
218a7a1495cSBarry Smith    TSComputeJacobian() is valid only for TS_NONLINEAR
219a7a1495cSBarry Smith 
220a7a1495cSBarry Smith    Level: developer
221a7a1495cSBarry Smith 
222a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
223a7a1495cSBarry Smith 
22494b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
225a7a1495cSBarry Smith @*/
22663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
227a7a1495cSBarry Smith {
228dfbe8321SBarry Smith   PetscErrorCode ierr;
229a7a1495cSBarry Smith 
230a7a1495cSBarry Smith   PetscFunctionBegin;
2314482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
2324482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
233c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
234a7a1495cSBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
23529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For TS_NONLINEAR only");
236a7a1495cSBarry Smith   }
237000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
238d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
239a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
240a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
241000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
242a7a1495cSBarry Smith     PetscStackPop;
243d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
244a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2454482741eSBarry Smith     PetscValidHeaderSpecific(*A,MAT_COOKIE,4);
2464482741eSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE,5);
247ef66eb69SBarry Smith   } else {
248ef66eb69SBarry Smith     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
249ef66eb69SBarry Smith     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
250ef66eb69SBarry Smith     if (*A != *B) {
251ef66eb69SBarry Smith       ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252ef66eb69SBarry Smith       ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
253ef66eb69SBarry Smith     }
254ef66eb69SBarry Smith   }
255a7a1495cSBarry Smith   PetscFunctionReturn(0);
256a7a1495cSBarry Smith }
257a7a1495cSBarry Smith 
2584a2ae208SSatish Balay #undef __FUNCT__
2594a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
260316643e7SJed Brown /*@
261d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
262d763cef2SBarry Smith 
263316643e7SJed Brown    Collective on TS and Vec
264316643e7SJed Brown 
265316643e7SJed Brown    Input Parameters:
266316643e7SJed Brown +  ts - the TS context
267316643e7SJed Brown .  t - current time
268316643e7SJed Brown -  x - state vector
269316643e7SJed Brown 
270316643e7SJed Brown    Output Parameter:
271316643e7SJed Brown .  y - right hand side
272316643e7SJed Brown 
273316643e7SJed Brown    Note:
274316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
275316643e7SJed Brown    is used internally within the nonlinear solvers.
276316643e7SJed Brown 
277316643e7SJed Brown    If the user did not provide a function but merely a matrix,
278d763cef2SBarry Smith    this routine applies the matrix.
279316643e7SJed Brown 
280316643e7SJed Brown    Level: developer
281316643e7SJed Brown 
282316643e7SJed Brown .keywords: TS, compute
283316643e7SJed Brown 
284316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction()
285316643e7SJed Brown @*/
286dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
287d763cef2SBarry Smith {
288dfbe8321SBarry Smith   PetscErrorCode ierr;
289d763cef2SBarry Smith 
290d763cef2SBarry Smith   PetscFunctionBegin;
2914482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
292316643e7SJed Brown   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
293316643e7SJed Brown   PetscValidHeaderSpecific(y,VEC_COOKIE,4);
294d763cef2SBarry Smith 
295d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
296000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
297d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
298000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
299d763cef2SBarry Smith     PetscStackPop;
3007c922b88SBarry Smith   } else {
301000e7ae3SMatthew Knepley     if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
302d763cef2SBarry Smith       MatStructure flg;
303d763cef2SBarry Smith       PetscStackPush("TS user right-hand-side matrix function");
3048beb423aSHong Zhang       ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
305d763cef2SBarry Smith       PetscStackPop;
306d763cef2SBarry Smith     }
3078beb423aSHong Zhang     ierr = MatMult(ts->Arhs,x,y);CHKERRQ(ierr);
3087c922b88SBarry Smith   }
309d763cef2SBarry Smith 
310d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
311d763cef2SBarry Smith 
312d763cef2SBarry Smith   PetscFunctionReturn(0);
313d763cef2SBarry Smith }
314d763cef2SBarry Smith 
3154a2ae208SSatish Balay #undef __FUNCT__
316316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction"
317316643e7SJed Brown /*@
318316643e7SJed Brown    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
319316643e7SJed Brown 
320316643e7SJed Brown    Collective on TS and Vec
321316643e7SJed Brown 
322316643e7SJed Brown    Input Parameters:
323316643e7SJed Brown +  ts - the TS context
324316643e7SJed Brown .  t - current time
325316643e7SJed Brown .  X - state vector
326316643e7SJed Brown -  Xdot - time derivative of state vector
327316643e7SJed Brown 
328316643e7SJed Brown    Output Parameter:
329316643e7SJed Brown .  Y - right hand side
330316643e7SJed Brown 
331316643e7SJed Brown    Note:
332316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
333316643e7SJed Brown    is used internally within the nonlinear solvers.
334316643e7SJed Brown 
335316643e7SJed Brown    If the user did did not write their equations in implicit form, this
336316643e7SJed Brown    function recasts them in implicit form.
337316643e7SJed Brown 
338316643e7SJed Brown    Level: developer
339316643e7SJed Brown 
340316643e7SJed Brown .keywords: TS, compute
341316643e7SJed Brown 
342316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction()
343316643e7SJed Brown @*/
344316643e7SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y)
345316643e7SJed Brown {
346316643e7SJed Brown   PetscErrorCode ierr;
347316643e7SJed Brown 
348316643e7SJed Brown   PetscFunctionBegin;
349316643e7SJed Brown   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
350316643e7SJed Brown   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
351316643e7SJed Brown   PetscValidHeaderSpecific(Xdot,VEC_COOKIE,4);
352316643e7SJed Brown   PetscValidHeaderSpecific(Y,VEC_COOKIE,5);
353316643e7SJed Brown 
354316643e7SJed Brown   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
355316643e7SJed Brown   if (ts->ops->ifunction) {
356316643e7SJed Brown     PetscStackPush("TS user implicit function");
357316643e7SJed Brown     ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr);
358316643e7SJed Brown     PetscStackPop;
359316643e7SJed Brown   } else {
360316643e7SJed Brown     PetscScalar *y;
361316643e7SJed Brown     PetscInt     i,n;
362316643e7SJed Brown 
363316643e7SJed Brown     if (ts->ops->rhsfunction) {
364316643e7SJed Brown       PetscStackPush("TS user right-hand-side function");
365316643e7SJed Brown       ierr = (*ts->ops->rhsfunction)(ts,t,X,Y,ts->funP);CHKERRQ(ierr);
366316643e7SJed Brown       PetscStackPop;
367316643e7SJed Brown     } else {
368316643e7SJed Brown       if (ts->ops->rhsmatrix) { /* assemble matrix for this timestep */
369316643e7SJed Brown         MatStructure flg;
370316643e7SJed Brown         PetscStackPush("TS user right-hand-side matrix function");
371316643e7SJed Brown         ierr = (*ts->ops->rhsmatrix)(ts,t,&ts->Arhs,&ts->B,&flg,ts->jacP);CHKERRQ(ierr);
372316643e7SJed Brown         PetscStackPop;
373316643e7SJed Brown       }
374316643e7SJed Brown       ierr = MatMult(ts->Arhs,X,Y);CHKERRQ(ierr);
375316643e7SJed Brown     }
376316643e7SJed Brown 
377316643e7SJed Brown     /* Convert to implicit form */
378316643e7SJed Brown     ierr = VecGetLocalSize(Y,&n);CHKERRQ(ierr);
379316643e7SJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
380316643e7SJed Brown     for (i=0; i<n; i++) y[i] = 1.0 - y[i];
381316643e7SJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
382316643e7SJed Brown   }
383316643e7SJed Brown   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
384316643e7SJed Brown   PetscFunctionReturn(0);
385316643e7SJed Brown }
386316643e7SJed Brown 
387316643e7SJed Brown #undef __FUNCT__
388316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian"
389316643e7SJed Brown /*@
390316643e7SJed Brown    TSComputeIJacobian - Evaluates the Jacobian of the DAE
391316643e7SJed Brown 
392316643e7SJed Brown    Collective on TS and Vec
393316643e7SJed Brown 
394316643e7SJed Brown    Input
395316643e7SJed Brown       Input Parameters:
396316643e7SJed Brown +  ts - the TS context
397316643e7SJed Brown .  t - current timestep
398316643e7SJed Brown .  X - state vector
399316643e7SJed Brown .  Xdot - time derivative of state vector
400316643e7SJed Brown -  shift - shift to apply, see note below
401316643e7SJed Brown 
402316643e7SJed Brown    Output Parameters:
403316643e7SJed Brown +  A - Jacobian matrix
404316643e7SJed Brown .  B - optional preconditioning matrix
405316643e7SJed Brown -  flag - flag indicating matrix structure
406316643e7SJed Brown 
407316643e7SJed Brown    Notes:
408316643e7SJed Brown    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
409316643e7SJed Brown 
410316643e7SJed Brown    dF/dX + shift*dF/dXdot
411316643e7SJed Brown 
412316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
413316643e7SJed Brown    is used internally within the nonlinear solvers.
414316643e7SJed Brown 
415316643e7SJed Brown    TSComputeIJacobian() is valid only for TS_NONLINEAR
416316643e7SJed Brown 
417316643e7SJed Brown    Level: developer
418316643e7SJed Brown 
419316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix
420316643e7SJed Brown 
421316643e7SJed Brown .seealso:  TSSetIJacobian()
422*63495f91SJed Brown @*/
423316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg)
424316643e7SJed Brown {
425316643e7SJed Brown   PetscErrorCode ierr;
426316643e7SJed Brown 
427316643e7SJed Brown   PetscFunctionBegin;
428316643e7SJed Brown   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
429316643e7SJed Brown   PetscValidHeaderSpecific(X,VEC_COOKIE,3);
430316643e7SJed Brown   PetscValidHeaderSpecific(Xdot,VEC_COOKIE,4);
431316643e7SJed Brown   PetscValidPointer(A,6);
432316643e7SJed Brown   PetscValidHeaderSpecific(*A,MAT_COOKIE,6);
433316643e7SJed Brown   PetscValidPointer(B,7);
434316643e7SJed Brown   PetscValidHeaderSpecific(*B,MAT_COOKIE,7);
435316643e7SJed Brown   PetscValidPointer(flg,8);
436316643e7SJed Brown 
437316643e7SJed Brown   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
438316643e7SJed Brown   if (ts->ops->ijacobian) {
439316643e7SJed Brown     PetscStackPush("TS user implicit Jacobian");
440316643e7SJed Brown     ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr);
441316643e7SJed Brown     PetscStackPop;
442316643e7SJed Brown   } else {
443316643e7SJed Brown     if (ts->ops->rhsjacobian) {
444316643e7SJed Brown       PetscStackPush("TS user right-hand-side Jacobian");
445316643e7SJed Brown       ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
446316643e7SJed Brown       PetscStackPop;
447316643e7SJed Brown     } else {
448316643e7SJed Brown       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
449316643e7SJed Brown       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
450316643e7SJed Brown       if (*A != *B) {
451316643e7SJed Brown         ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
452316643e7SJed Brown         ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
453316643e7SJed Brown       }
454316643e7SJed Brown     }
455316643e7SJed Brown 
456316643e7SJed Brown     /* Convert to implicit form */
457316643e7SJed Brown     /* inefficient because these operations will normally traverse all matrix elements twice */
458316643e7SJed Brown     ierr = MatScale(*A,-1);CHKERRQ(ierr);
459316643e7SJed Brown     ierr = MatShift(*A,shift);CHKERRQ(ierr);
460316643e7SJed Brown     if (*A != *B) {
461316643e7SJed Brown       ierr = MatScale(*B,-1);CHKERRQ(ierr);
462316643e7SJed Brown       ierr = MatShift(*B,shift);CHKERRQ(ierr);
463316643e7SJed Brown     }
464316643e7SJed Brown   }
465316643e7SJed Brown   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
466316643e7SJed Brown   PetscFunctionReturn(0);
467316643e7SJed Brown }
468316643e7SJed Brown 
469316643e7SJed Brown #undef __FUNCT__
4704a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
471d763cef2SBarry Smith /*@C
472d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
473d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
474d763cef2SBarry Smith 
475d763cef2SBarry Smith     Collective on TS
476d763cef2SBarry Smith 
477d763cef2SBarry Smith     Input Parameters:
478d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
479d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
480d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
481d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
482d763cef2SBarry Smith 
483d763cef2SBarry Smith     Calling sequence of func:
48487828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
485d763cef2SBarry Smith 
486d763cef2SBarry Smith +   t - current timestep
487d763cef2SBarry Smith .   u - input vector
488d763cef2SBarry Smith .   F - function vector
489d763cef2SBarry Smith -   ctx - [optional] user-defined function context
490d763cef2SBarry Smith 
491d763cef2SBarry Smith     Important:
49276f2fa84SHong Zhang     The user MUST call either this routine or TSSetMatrices().
493d763cef2SBarry Smith 
494d763cef2SBarry Smith     Level: beginner
495d763cef2SBarry Smith 
496d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
497d763cef2SBarry Smith 
49876f2fa84SHong Zhang .seealso: TSSetMatrices()
499d763cef2SBarry Smith @*/
50063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
501d763cef2SBarry Smith {
502d763cef2SBarry Smith   PetscFunctionBegin;
503d763cef2SBarry Smith 
5044482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
505d763cef2SBarry Smith   if (ts->problem_type == TS_LINEAR) {
50629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot set function for linear problem");
507d763cef2SBarry Smith   }
508000e7ae3SMatthew Knepley   ts->ops->rhsfunction = f;
509d763cef2SBarry Smith   ts->funP             = ctx;
510d763cef2SBarry Smith   PetscFunctionReturn(0);
511d763cef2SBarry Smith }
512d763cef2SBarry Smith 
5134a2ae208SSatish Balay #undef __FUNCT__
51495f0b562SHong Zhang #define __FUNCT__ "TSSetMatrices"
51595f0b562SHong Zhang /*@C
51695f0b562SHong Zhang    TSSetMatrices - Sets the functions to compute the matrices Alhs and Arhs,
51795f0b562SHong Zhang    where Alhs(t) U_t = Arhs(t) U.
51895f0b562SHong Zhang 
51995f0b562SHong Zhang    Collective on TS
52095f0b562SHong Zhang 
52195f0b562SHong Zhang    Input Parameters:
52295f0b562SHong Zhang +  ts   - the TS context obtained from TSCreate()
52395f0b562SHong Zhang .  Arhs - matrix
52495f0b562SHong Zhang .  frhs - the matrix evaluation routine for Arhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
52595f0b562SHong Zhang           if Arhs is not a function of t.
52695f0b562SHong Zhang .  Alhs - matrix or PETSC_NULL if Alhs is an indentity matrix.
52795f0b562SHong Zhang .  flhs - the matrix evaluation routine for Alhs; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran)
52895f0b562SHong Zhang           if Alhs is not a function of t.
52995f0b562SHong Zhang .  flag - flag indicating information about the matrix structure of Arhs and Alhs.
53095f0b562SHong Zhang           The available options are
53195f0b562SHong Zhang             SAME_NONZERO_PATTERN - Alhs has the same nonzero structure as Arhs
53295f0b562SHong Zhang             DIFFERENT_NONZERO_PATTERN - Alhs has different nonzero structure as Arhs
53395f0b562SHong Zhang -  ctx  - [optional] user-defined context for private data for the
53495f0b562SHong Zhang           matrix evaluation routine (may be PETSC_NULL)
53595f0b562SHong Zhang 
53695f0b562SHong Zhang    Calling sequence of func:
53795f0b562SHong Zhang $     func(TS ts,PetscReal t,Mat *A,Mat *B,PetscInt *flag,void *ctx);
53895f0b562SHong Zhang 
53995f0b562SHong Zhang +  t - current timestep
54095f0b562SHong Zhang .  A - matrix A, where U_t = A(t) U
54195f0b562SHong Zhang .  B - preconditioner matrix, usually the same as A
54295f0b562SHong Zhang .  flag - flag indicating information about the preconditioner matrix
54395f0b562SHong Zhang           structure (same as flag in KSPSetOperators())
54495f0b562SHong Zhang -  ctx - [optional] user-defined context for matrix evaluation routine
54595f0b562SHong Zhang 
54695f0b562SHong Zhang    Notes:
54795f0b562SHong Zhang    The routine func() takes Mat* as the matrix arguments rather than Mat.
54895f0b562SHong Zhang    This allows the matrix evaluation routine to replace Arhs or Alhs with a
54995f0b562SHong Zhang    completely new new matrix structure (not just different matrix elements)
55095f0b562SHong Zhang    when appropriate, for instance, if the nonzero structure is changing
55195f0b562SHong Zhang    throughout the global iterations.
55295f0b562SHong Zhang 
55395f0b562SHong Zhang    Important:
55495f0b562SHong Zhang    The user MUST call either this routine or TSSetRHSFunction().
55595f0b562SHong Zhang 
55695f0b562SHong Zhang    Level: beginner
55795f0b562SHong Zhang 
55895f0b562SHong Zhang .keywords: TS, timestep, set, matrix
55995f0b562SHong Zhang 
56095f0b562SHong Zhang .seealso: TSSetRHSFunction()
56195f0b562SHong Zhang @*/
56295f0b562SHong 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)
56395f0b562SHong Zhang {
56495f0b562SHong Zhang   PetscFunctionBegin;
56595f0b562SHong Zhang   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
56692af4f6aSHong Zhang   if (Arhs){
56795f0b562SHong Zhang     PetscValidHeaderSpecific(Arhs,MAT_COOKIE,2);
56895f0b562SHong Zhang     PetscCheckSameComm(ts,1,Arhs,2);
56995f0b562SHong Zhang     ts->Arhs           = Arhs;
57092af4f6aSHong Zhang     ts->ops->rhsmatrix = frhs;
57192af4f6aSHong Zhang   }
57292af4f6aSHong Zhang   if (Alhs){
57392af4f6aSHong Zhang     PetscValidHeaderSpecific(Alhs,MAT_COOKIE,4);
57492af4f6aSHong Zhang     PetscCheckSameComm(ts,1,Alhs,4);
57595f0b562SHong Zhang     ts->Alhs           = Alhs;
57692af4f6aSHong Zhang     ts->ops->lhsmatrix = flhs;
57792af4f6aSHong Zhang   }
57892af4f6aSHong Zhang 
57992af4f6aSHong Zhang   ts->jacP           = ctx;
58095f0b562SHong Zhang   ts->matflg         = flag;
58195f0b562SHong Zhang   PetscFunctionReturn(0);
58295f0b562SHong Zhang }
583d763cef2SBarry Smith 
584aa644b49SHong Zhang #undef __FUNCT__
585cda39b92SHong Zhang #define __FUNCT__ "TSGetMatrices"
586cda39b92SHong Zhang /*@C
587cda39b92SHong Zhang    TSGetMatrices - Returns the matrices Arhs and Alhs at the present timestep,
588cda39b92SHong Zhang    where Alhs(t) U_t = Arhs(t) U.
589cda39b92SHong Zhang 
590cda39b92SHong Zhang    Not Collective, but parallel objects are returned if TS is parallel
591cda39b92SHong Zhang 
592cda39b92SHong Zhang    Input Parameter:
593cda39b92SHong Zhang .  ts  - The TS context obtained from TSCreate()
594cda39b92SHong Zhang 
595cda39b92SHong Zhang    Output Parameters:
596cda39b92SHong Zhang +  Arhs - The right-hand side matrix
597cda39b92SHong Zhang .  Alhs - The left-hand side matrix
598cda39b92SHong Zhang -  ctx - User-defined context for matrix evaluation routine
599cda39b92SHong Zhang 
600cda39b92SHong Zhang    Notes: You can pass in PETSC_NULL for any return argument you do not need.
601cda39b92SHong Zhang 
602cda39b92SHong Zhang    Level: intermediate
603cda39b92SHong Zhang 
604cda39b92SHong Zhang .seealso: TSSetMatrices(), TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian()
605cda39b92SHong Zhang 
606cda39b92SHong Zhang .keywords: TS, timestep, get, matrix
607cda39b92SHong Zhang 
608cda39b92SHong Zhang @*/
609cda39b92SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSGetMatrices(TS ts,Mat *Arhs,Mat *Alhs,void **ctx)
610cda39b92SHong Zhang {
611cda39b92SHong Zhang   PetscFunctionBegin;
612cda39b92SHong Zhang   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
613cda39b92SHong Zhang   if (Arhs) *Arhs = ts->Arhs;
614cda39b92SHong Zhang   if (Alhs) *Alhs = ts->Alhs;
615cda39b92SHong Zhang   if (ctx)  *ctx = ts->jacP;
616cda39b92SHong Zhang   PetscFunctionReturn(0);
617cda39b92SHong Zhang }
618cda39b92SHong Zhang 
619cda39b92SHong Zhang #undef __FUNCT__
6204a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
621d763cef2SBarry Smith /*@C
622d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
623d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
62476f2fa84SHong Zhang    Use TSSetMatrices() for linear problems.
625d763cef2SBarry Smith 
626d763cef2SBarry Smith    Collective on TS
627d763cef2SBarry Smith 
628d763cef2SBarry Smith    Input Parameters:
629d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
630d763cef2SBarry Smith .  A   - Jacobian matrix
631d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
632d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
633d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
634d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
635d763cef2SBarry Smith 
636d763cef2SBarry Smith    Calling sequence of func:
63787828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
638d763cef2SBarry Smith 
639d763cef2SBarry Smith +  t - current timestep
640d763cef2SBarry Smith .  u - input vector
641d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
642d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
643d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
64494b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
645d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
646d763cef2SBarry Smith 
647d763cef2SBarry Smith    Notes:
64894b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
649d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
650d763cef2SBarry Smith 
651d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
652d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
65356335db2SHong Zhang    completely new matrix structure (not just different matrix elements)
654d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
655d763cef2SBarry Smith    throughout the global iterations.
656d763cef2SBarry Smith 
657d763cef2SBarry Smith    Level: beginner
658d763cef2SBarry Smith 
659d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
660d763cef2SBarry Smith 
661d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
66276f2fa84SHong Zhang           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices()
663d763cef2SBarry Smith 
664d763cef2SBarry Smith @*/
66563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetRHSJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
666d763cef2SBarry Smith {
667d763cef2SBarry Smith   PetscFunctionBegin;
6684482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
6694482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
6704482741eSBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE,3);
671c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,A,2);
672c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,B,3);
673d763cef2SBarry Smith   if (ts->problem_type != TS_NONLINEAR) {
67476f2fa84SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONG,"Not for linear problems; use TSSetMatrices()");
675d763cef2SBarry Smith   }
676d763cef2SBarry Smith 
677000e7ae3SMatthew Knepley   ts->ops->rhsjacobian = f;
678d763cef2SBarry Smith   ts->jacP             = ctx;
6798beb423aSHong Zhang   ts->Arhs             = A;
680d763cef2SBarry Smith   ts->B                = B;
681d763cef2SBarry Smith   PetscFunctionReturn(0);
682d763cef2SBarry Smith }
683d763cef2SBarry Smith 
684316643e7SJed Brown 
685316643e7SJed Brown #undef __FUNCT__
686316643e7SJed Brown #define __FUNCT__ "TSSetIFunction"
687316643e7SJed Brown /*@C
688316643e7SJed Brown    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
689316643e7SJed Brown 
690316643e7SJed Brown    Collective on TS
691316643e7SJed Brown 
692316643e7SJed Brown    Input Parameters:
693316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
694316643e7SJed Brown .  f   - the function evaluation routine
695316643e7SJed Brown -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
696316643e7SJed Brown 
697316643e7SJed Brown    Calling sequence of f:
698316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
699316643e7SJed Brown 
700316643e7SJed Brown +  t   - time at step/stage being solved
701316643e7SJed Brown .  u   - state vector
702316643e7SJed Brown .  u_t - time derivative of state vector
703316643e7SJed Brown .  F   - function vector
704316643e7SJed Brown -  ctx - [optional] user-defined context for matrix evaluation routine
705316643e7SJed Brown 
706316643e7SJed Brown    Important:
707316643e7SJed Brown    The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices().  This routine must be used when not solving an ODE.
708316643e7SJed Brown 
709316643e7SJed Brown    Level: beginner
710316643e7SJed Brown 
711316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian
712316643e7SJed Brown 
713316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian()
714316643e7SJed Brown @*/
715316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIFunction(TS ts,TSIFunction f,void *ctx)
716316643e7SJed Brown {
717316643e7SJed Brown 
718316643e7SJed Brown   PetscFunctionBegin;
719316643e7SJed Brown   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
720316643e7SJed Brown   ts->ops->ifunction = f;
721316643e7SJed Brown   ts->funP           = ctx;
722316643e7SJed Brown   PetscFunctionReturn(0);
723316643e7SJed Brown }
724316643e7SJed Brown 
725316643e7SJed Brown #undef __FUNCT__
726316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian"
727316643e7SJed Brown /*@C
728316643e7SJed Brown    TSSetIJacobian - Set the function to compute the Jacobian of
729316643e7SJed Brown    G(U) = F(t,U,U0+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
730316643e7SJed Brown 
731316643e7SJed Brown    Collective on TS
732316643e7SJed Brown 
733316643e7SJed Brown    Input Parameters:
734316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
735316643e7SJed Brown .  A   - Jacobian matrix
736316643e7SJed Brown .  B   - preconditioning matrix for A (may be same as A)
737316643e7SJed Brown .  f   - the Jacobian evaluation routine
738316643e7SJed Brown -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
739316643e7SJed Brown 
740316643e7SJed Brown    Calling sequence of f:
741316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
742316643e7SJed Brown 
743316643e7SJed Brown +  t    - time at step/stage being solved
744316643e7SJed Brown .  u    - state vector
745316643e7SJed Brown .  u_t  - time derivative of state vector
746316643e7SJed Brown .  a    - shift
747316643e7SJed Brown .  A    - Jacobian of G(U) = F(t,U,U0+a*U), equivalent to dF/dU + a*dF/dU_t
748316643e7SJed Brown .  B    - preconditioning matrix for A, may be same as A
749316643e7SJed Brown .  flag - flag indicating information about the preconditioner matrix
750316643e7SJed Brown           structure (same as flag in KSPSetOperators())
751316643e7SJed Brown -  ctx  - [optional] user-defined context for matrix evaluation routine
752316643e7SJed Brown 
753316643e7SJed Brown    Notes:
754316643e7SJed Brown    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
755316643e7SJed Brown 
756316643e7SJed Brown    Level: beginner
757316643e7SJed Brown 
758316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian
759316643e7SJed Brown 
760316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian()
761316643e7SJed Brown 
762316643e7SJed Brown @*/
763316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
764316643e7SJed Brown {
765316643e7SJed Brown   PetscErrorCode ierr;
766316643e7SJed Brown 
767316643e7SJed Brown   PetscFunctionBegin;
768316643e7SJed Brown   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
769316643e7SJed Brown   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
770316643e7SJed Brown   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
771316643e7SJed Brown   if (A) PetscCheckSameComm(ts,1,A,2);
772316643e7SJed Brown   if (B) PetscCheckSameComm(ts,1,B,3);
773316643e7SJed Brown   if (f)   ts->ops->ijacobian = f;
774316643e7SJed Brown   if (ctx) ts->jacP             = ctx;
775316643e7SJed Brown   if (A) {
776316643e7SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
777316643e7SJed Brown     if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr);}
778316643e7SJed Brown     ts->A = A;
779316643e7SJed Brown   }
780316643e7SJed Brown   if (B) {
781316643e7SJed Brown     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
782316643e7SJed Brown     if (ts->B) {ierr = MatDestroy(ts->B);CHKERRQ(ierr);}
783316643e7SJed Brown     ts->B = B;
784316643e7SJed Brown   }
785316643e7SJed Brown   PetscFunctionReturn(0);
786316643e7SJed Brown }
787316643e7SJed Brown 
7884a2ae208SSatish Balay #undef __FUNCT__
7894a2ae208SSatish Balay #define __FUNCT__ "TSView"
7907e2c5f70SBarry Smith /*@C
791d763cef2SBarry Smith     TSView - Prints the TS data structure.
792d763cef2SBarry Smith 
7934c49b128SBarry Smith     Collective on TS
794d763cef2SBarry Smith 
795d763cef2SBarry Smith     Input Parameters:
796d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
797d763cef2SBarry Smith -   viewer - visualization context
798d763cef2SBarry Smith 
799d763cef2SBarry Smith     Options Database Key:
800d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
801d763cef2SBarry Smith 
802d763cef2SBarry Smith     Notes:
803d763cef2SBarry Smith     The available visualization contexts include
804b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
805b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
806d763cef2SBarry Smith          output where only the first processor opens
807d763cef2SBarry Smith          the file.  All other processors send their
808d763cef2SBarry Smith          data to the first processor to print.
809d763cef2SBarry Smith 
810d763cef2SBarry Smith     The user can open an alternative visualization context with
811b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
812d763cef2SBarry Smith 
813d763cef2SBarry Smith     Level: beginner
814d763cef2SBarry Smith 
815d763cef2SBarry Smith .keywords: TS, timestep, view
816d763cef2SBarry Smith 
817b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
818d763cef2SBarry Smith @*/
81963dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSView(TS ts,PetscViewer viewer)
820d763cef2SBarry Smith {
821dfbe8321SBarry Smith   PetscErrorCode ierr;
822a313700dSBarry Smith   const TSType   type;
82332077d6dSBarry Smith   PetscTruth     iascii,isstring;
824d763cef2SBarry Smith 
825d763cef2SBarry Smith   PetscFunctionBegin;
8264482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
8273050cee2SBarry Smith   if (!viewer) {
8287adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
8293050cee2SBarry Smith   }
8304482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
831c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
832fd16b177SBarry Smith 
83332077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
834b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
83532077d6dSBarry Smith   if (iascii) {
836b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr);
837a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
838454a90a3SBarry Smith     if (type) {
839b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
840184914b5SBarry Smith     } else {
841b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
842184914b5SBarry Smith     }
843000e7ae3SMatthew Knepley     if (ts->ops->view) {
844b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
845000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
846b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
847d763cef2SBarry Smith     }
84877431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
849a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
850d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
85177431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
852d763cef2SBarry Smith     }
85377431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
8540f5bd95cSBarry Smith   } else if (isstring) {
855a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
856b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
857d763cef2SBarry Smith   }
858b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
85994b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPView(ts->ksp,viewer);CHKERRQ(ierr);}
860d763cef2SBarry Smith   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
861b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
862d763cef2SBarry Smith   PetscFunctionReturn(0);
863d763cef2SBarry Smith }
864d763cef2SBarry Smith 
865d763cef2SBarry Smith 
8664a2ae208SSatish Balay #undef __FUNCT__
8674a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
868d763cef2SBarry Smith /*@C
869d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
870d763cef2SBarry Smith    the timesteppers.
871d763cef2SBarry Smith 
872d763cef2SBarry Smith    Collective on TS
873d763cef2SBarry Smith 
874d763cef2SBarry Smith    Input Parameters:
875d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
876d763cef2SBarry Smith -  usrP - optional user context
877d763cef2SBarry Smith 
878d763cef2SBarry Smith    Level: intermediate
879d763cef2SBarry Smith 
880d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
881d763cef2SBarry Smith 
882d763cef2SBarry Smith .seealso: TSGetApplicationContext()
883d763cef2SBarry Smith @*/
88463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetApplicationContext(TS ts,void *usrP)
885d763cef2SBarry Smith {
886d763cef2SBarry Smith   PetscFunctionBegin;
8874482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
888d763cef2SBarry Smith   ts->user = usrP;
889d763cef2SBarry Smith   PetscFunctionReturn(0);
890d763cef2SBarry Smith }
891d763cef2SBarry Smith 
8924a2ae208SSatish Balay #undef __FUNCT__
8934a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
894d763cef2SBarry Smith /*@C
895d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
896d763cef2SBarry Smith     timestepper.
897d763cef2SBarry Smith 
898d763cef2SBarry Smith     Not Collective
899d763cef2SBarry Smith 
900d763cef2SBarry Smith     Input Parameter:
901d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
902d763cef2SBarry Smith 
903d763cef2SBarry Smith     Output Parameter:
904d763cef2SBarry Smith .   usrP - user context
905d763cef2SBarry Smith 
906d763cef2SBarry Smith     Level: intermediate
907d763cef2SBarry Smith 
908d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
909d763cef2SBarry Smith 
910d763cef2SBarry Smith .seealso: TSSetApplicationContext()
911d763cef2SBarry Smith @*/
91263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetApplicationContext(TS ts,void **usrP)
913d763cef2SBarry Smith {
914d763cef2SBarry Smith   PetscFunctionBegin;
9154482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
916d763cef2SBarry Smith   *usrP = ts->user;
917d763cef2SBarry Smith   PetscFunctionReturn(0);
918d763cef2SBarry Smith }
919d763cef2SBarry Smith 
9204a2ae208SSatish Balay #undef __FUNCT__
9214a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
922d763cef2SBarry Smith /*@
923d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
924d763cef2SBarry Smith 
925d763cef2SBarry Smith    Not Collective
926d763cef2SBarry Smith 
927d763cef2SBarry Smith    Input Parameter:
928d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
929d763cef2SBarry Smith 
930d763cef2SBarry Smith    Output Parameter:
931d763cef2SBarry Smith .  iter - number steps so far
932d763cef2SBarry Smith 
933d763cef2SBarry Smith    Level: intermediate
934d763cef2SBarry Smith 
935d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
936d763cef2SBarry Smith @*/
93763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStepNumber(TS ts,PetscInt* iter)
938d763cef2SBarry Smith {
939d763cef2SBarry Smith   PetscFunctionBegin;
9404482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
9414482741eSBarry Smith   PetscValidIntPointer(iter,2);
942d763cef2SBarry Smith   *iter = ts->steps;
943d763cef2SBarry Smith   PetscFunctionReturn(0);
944d763cef2SBarry Smith }
945d763cef2SBarry Smith 
9464a2ae208SSatish Balay #undef __FUNCT__
9474a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
948d763cef2SBarry Smith /*@
949d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
950d763cef2SBarry Smith    as well as the initial time.
951d763cef2SBarry Smith 
952d763cef2SBarry Smith    Collective on TS
953d763cef2SBarry Smith 
954d763cef2SBarry Smith    Input Parameters:
955d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
956d763cef2SBarry Smith .  initial_time - the initial time
957d763cef2SBarry Smith -  time_step - the size of the timestep
958d763cef2SBarry Smith 
959d763cef2SBarry Smith    Level: intermediate
960d763cef2SBarry Smith 
961d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
962d763cef2SBarry Smith 
963d763cef2SBarry Smith .keywords: TS, set, initial, timestep
964d763cef2SBarry Smith @*/
96563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
966d763cef2SBarry Smith {
967d763cef2SBarry Smith   PetscFunctionBegin;
9684482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
969d763cef2SBarry Smith   ts->time_step         = time_step;
970d763cef2SBarry Smith   ts->initial_time_step = time_step;
971d763cef2SBarry Smith   ts->ptime             = initial_time;
972d763cef2SBarry Smith   PetscFunctionReturn(0);
973d763cef2SBarry Smith }
974d763cef2SBarry Smith 
9754a2ae208SSatish Balay #undef __FUNCT__
9764a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
977d763cef2SBarry Smith /*@
978d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
979d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
980d763cef2SBarry Smith 
981d763cef2SBarry Smith    Collective on TS
982d763cef2SBarry Smith 
983d763cef2SBarry Smith    Input Parameters:
984d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
985d763cef2SBarry Smith -  time_step - the size of the timestep
986d763cef2SBarry Smith 
987d763cef2SBarry Smith    Level: intermediate
988d763cef2SBarry Smith 
989d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
990d763cef2SBarry Smith 
991d763cef2SBarry Smith .keywords: TS, set, timestep
992d763cef2SBarry Smith @*/
99363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetTimeStep(TS ts,PetscReal time_step)
994d763cef2SBarry Smith {
995d763cef2SBarry Smith   PetscFunctionBegin;
9964482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
997d763cef2SBarry Smith   ts->time_step = time_step;
998d763cef2SBarry Smith   PetscFunctionReturn(0);
999d763cef2SBarry Smith }
1000d763cef2SBarry Smith 
10014a2ae208SSatish Balay #undef __FUNCT__
10024a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
1003d763cef2SBarry Smith /*@
1004d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
1005d763cef2SBarry Smith 
1006d763cef2SBarry Smith    Not Collective
1007d763cef2SBarry Smith 
1008d763cef2SBarry Smith    Input Parameter:
1009d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1010d763cef2SBarry Smith 
1011d763cef2SBarry Smith    Output Parameter:
1012d763cef2SBarry Smith .  dt - the current timestep size
1013d763cef2SBarry Smith 
1014d763cef2SBarry Smith    Level: intermediate
1015d763cef2SBarry Smith 
1016d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1017d763cef2SBarry Smith 
1018d763cef2SBarry Smith .keywords: TS, get, timestep
1019d763cef2SBarry Smith @*/
102063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTimeStep(TS ts,PetscReal* dt)
1021d763cef2SBarry Smith {
1022d763cef2SBarry Smith   PetscFunctionBegin;
10234482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
10244482741eSBarry Smith   PetscValidDoublePointer(dt,2);
1025d763cef2SBarry Smith   *dt = ts->time_step;
1026d763cef2SBarry Smith   PetscFunctionReturn(0);
1027d763cef2SBarry Smith }
1028d763cef2SBarry Smith 
10294a2ae208SSatish Balay #undef __FUNCT__
10304a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
1031d8e5e3e6SSatish Balay /*@
1032d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
1033d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
1034d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
1035d763cef2SBarry Smith    the solution at the next timestep has been calculated.
1036d763cef2SBarry Smith 
1037d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
1038d763cef2SBarry Smith 
1039d763cef2SBarry Smith    Input Parameter:
1040d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1041d763cef2SBarry Smith 
1042d763cef2SBarry Smith    Output Parameter:
1043d763cef2SBarry Smith .  v - the vector containing the solution
1044d763cef2SBarry Smith 
1045d763cef2SBarry Smith    Level: intermediate
1046d763cef2SBarry Smith 
1047d763cef2SBarry Smith .seealso: TSGetTimeStep()
1048d763cef2SBarry Smith 
1049d763cef2SBarry Smith .keywords: TS, timestep, get, solution
1050d763cef2SBarry Smith @*/
105163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSolution(TS ts,Vec *v)
1052d763cef2SBarry Smith {
1053d763cef2SBarry Smith   PetscFunctionBegin;
10544482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
10554482741eSBarry Smith   PetscValidPointer(v,2);
1056d763cef2SBarry Smith   *v = ts->vec_sol_always;
1057d763cef2SBarry Smith   PetscFunctionReturn(0);
1058d763cef2SBarry Smith }
1059d763cef2SBarry Smith 
1060bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
10614a2ae208SSatish Balay #undef __FUNCT__
1062bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
1063d8e5e3e6SSatish Balay /*@
1064bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
1065d763cef2SBarry Smith 
1066bdad233fSMatthew Knepley   Not collective
1067d763cef2SBarry Smith 
1068bdad233fSMatthew Knepley   Input Parameters:
1069bdad233fSMatthew Knepley + ts   - The TS
1070bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1071d763cef2SBarry Smith .vb
1072d763cef2SBarry Smith          U_t = A U
1073d763cef2SBarry Smith          U_t = A(t) U
1074d763cef2SBarry Smith          U_t = F(t,U)
1075d763cef2SBarry Smith .ve
1076d763cef2SBarry Smith 
1077d763cef2SBarry Smith    Level: beginner
1078d763cef2SBarry Smith 
1079bdad233fSMatthew Knepley .keywords: TS, problem type
1080bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1081d763cef2SBarry Smith @*/
108263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetProblemType(TS ts, TSProblemType type)
1083a7cc72afSBarry Smith {
1084d763cef2SBarry Smith   PetscFunctionBegin;
10854482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1086bdad233fSMatthew Knepley   ts->problem_type = type;
1087d763cef2SBarry Smith   PetscFunctionReturn(0);
1088d763cef2SBarry Smith }
1089d763cef2SBarry Smith 
1090bdad233fSMatthew Knepley #undef __FUNCT__
1091bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
1092bdad233fSMatthew Knepley /*@C
1093bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
1094bdad233fSMatthew Knepley 
1095bdad233fSMatthew Knepley   Not collective
1096bdad233fSMatthew Knepley 
1097bdad233fSMatthew Knepley   Input Parameter:
1098bdad233fSMatthew Knepley . ts   - The TS
1099bdad233fSMatthew Knepley 
1100bdad233fSMatthew Knepley   Output Parameter:
1101bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1102bdad233fSMatthew Knepley .vb
1103bdad233fSMatthew Knepley          U_t = A U
1104bdad233fSMatthew Knepley          U_t = A(t) U
1105bdad233fSMatthew Knepley          U_t = F(t,U)
1106bdad233fSMatthew Knepley .ve
1107bdad233fSMatthew Knepley 
1108bdad233fSMatthew Knepley    Level: beginner
1109bdad233fSMatthew Knepley 
1110bdad233fSMatthew Knepley .keywords: TS, problem type
1111bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1112bdad233fSMatthew Knepley @*/
111363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetProblemType(TS ts, TSProblemType *type)
1114a7cc72afSBarry Smith {
1115bdad233fSMatthew Knepley   PetscFunctionBegin;
11164482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
11174482741eSBarry Smith   PetscValidIntPointer(type,2);
1118bdad233fSMatthew Knepley   *type = ts->problem_type;
1119bdad233fSMatthew Knepley   PetscFunctionReturn(0);
1120bdad233fSMatthew Knepley }
1121d763cef2SBarry Smith 
11224a2ae208SSatish Balay #undef __FUNCT__
11234a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
1124d763cef2SBarry Smith /*@
1125d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
1126d763cef2SBarry Smith    of a timestepper.
1127d763cef2SBarry Smith 
1128d763cef2SBarry Smith    Collective on TS
1129d763cef2SBarry Smith 
1130d763cef2SBarry Smith    Input Parameter:
1131d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1132d763cef2SBarry Smith 
1133d763cef2SBarry Smith    Notes:
1134d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
1135d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
1136d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
1137d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
1138d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
1139d763cef2SBarry Smith 
1140d763cef2SBarry Smith    Level: advanced
1141d763cef2SBarry Smith 
1142d763cef2SBarry Smith .keywords: TS, timestep, setup
1143d763cef2SBarry Smith 
1144d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
1145d763cef2SBarry Smith @*/
114663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetUp(TS ts)
1147d763cef2SBarry Smith {
1148dfbe8321SBarry Smith   PetscErrorCode ierr;
1149d763cef2SBarry Smith 
1150d763cef2SBarry Smith   PetscFunctionBegin;
11514482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
115229bbc08cSBarry Smith   if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
11537adad957SLisandro Dalcin   if (!((PetscObject)ts)->type_name) {
1154d763cef2SBarry Smith     ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr);
1155d763cef2SBarry Smith   }
1156000e7ae3SMatthew Knepley   ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1157d763cef2SBarry Smith   ts->setupcalled = 1;
1158d763cef2SBarry Smith   PetscFunctionReturn(0);
1159d763cef2SBarry Smith }
1160d763cef2SBarry Smith 
11614a2ae208SSatish Balay #undef __FUNCT__
11624a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
1163d8e5e3e6SSatish Balay /*@
1164d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
1165d763cef2SBarry Smith    with TSCreate().
1166d763cef2SBarry Smith 
1167d763cef2SBarry Smith    Collective on TS
1168d763cef2SBarry Smith 
1169d763cef2SBarry Smith    Input Parameter:
1170d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1171d763cef2SBarry Smith 
1172d763cef2SBarry Smith    Level: beginner
1173d763cef2SBarry Smith 
1174d763cef2SBarry Smith .keywords: TS, timestepper, destroy
1175d763cef2SBarry Smith 
1176d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
1177d763cef2SBarry Smith @*/
117863dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDestroy(TS ts)
1179d763cef2SBarry Smith {
11806849ba73SBarry Smith   PetscErrorCode ierr;
1181d763cef2SBarry Smith 
1182d763cef2SBarry Smith   PetscFunctionBegin;
11834482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
11847adad957SLisandro Dalcin   if (--((PetscObject)ts)->refct > 0) PetscFunctionReturn(0);
1185d763cef2SBarry Smith 
1186be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
11870f5bd95cSBarry Smith   ierr = PetscObjectDepublish(ts);CHKERRQ(ierr);
11888beb423aSHong Zhang   if (ts->A) {ierr = MatDestroy(ts->A);CHKERRQ(ierr)}
118994b7f48cSBarry Smith   if (ts->ksp) {ierr = KSPDestroy(ts->ksp);CHKERRQ(ierr);}
1190d763cef2SBarry Smith   if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);}
11911e3347e8SBarry Smith   if (ts->ops->destroy) {ierr = (*(ts)->ops->destroy)(ts);CHKERRQ(ierr);}
1192a6570f20SBarry Smith   ierr = TSMonitorCancel(ts);CHKERRQ(ierr);
1193a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1194d763cef2SBarry Smith   PetscFunctionReturn(0);
1195d763cef2SBarry Smith }
1196d763cef2SBarry Smith 
11974a2ae208SSatish Balay #undef __FUNCT__
11984a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
1199d8e5e3e6SSatish Balay /*@
1200d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1201d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
1202d763cef2SBarry Smith 
1203d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
1204d763cef2SBarry Smith 
1205d763cef2SBarry Smith    Input Parameter:
1206d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1207d763cef2SBarry Smith 
1208d763cef2SBarry Smith    Output Parameter:
1209d763cef2SBarry Smith .  snes - the nonlinear solver context
1210d763cef2SBarry Smith 
1211d763cef2SBarry Smith    Notes:
1212d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
1213d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
121494b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
1215d763cef2SBarry Smith 
1216d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
1217d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
1218d763cef2SBarry Smith 
1219d763cef2SBarry Smith    Level: beginner
1220d763cef2SBarry Smith 
1221d763cef2SBarry Smith .keywords: timestep, get, SNES
1222d763cef2SBarry Smith @*/
122363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetSNES(TS ts,SNES *snes)
1224d763cef2SBarry Smith {
1225d763cef2SBarry Smith   PetscFunctionBegin;
12264482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
12274482741eSBarry Smith   PetscValidPointer(snes,2);
1228447600ffSHong Zhang   if (((PetscObject)ts)->type_name == PETSC_NULL)
1229447600ffSHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"SNES is not created yet. Call TSSetType() first");
123094b7f48cSBarry Smith   if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Nonlinear only; use TSGetKSP()");
1231d763cef2SBarry Smith   *snes = ts->snes;
1232d763cef2SBarry Smith   PetscFunctionReturn(0);
1233d763cef2SBarry Smith }
1234d763cef2SBarry Smith 
12354a2ae208SSatish Balay #undef __FUNCT__
123694b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
1237d8e5e3e6SSatish Balay /*@
123894b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
1239d763cef2SBarry Smith    a TS (timestepper) context.
1240d763cef2SBarry Smith 
124194b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1242d763cef2SBarry Smith 
1243d763cef2SBarry Smith    Input Parameter:
1244d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1245d763cef2SBarry Smith 
1246d763cef2SBarry Smith    Output Parameter:
124794b7f48cSBarry Smith .  ksp - the nonlinear solver context
1248d763cef2SBarry Smith 
1249d763cef2SBarry Smith    Notes:
125094b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1251d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1252d763cef2SBarry Smith    KSP and PC contexts as well.
1253d763cef2SBarry Smith 
125494b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
125594b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1256d763cef2SBarry Smith 
1257d763cef2SBarry Smith    Level: beginner
1258d763cef2SBarry Smith 
125994b7f48cSBarry Smith .keywords: timestep, get, KSP
1260d763cef2SBarry Smith @*/
126163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetKSP(TS ts,KSP *ksp)
1262d763cef2SBarry Smith {
1263d763cef2SBarry Smith   PetscFunctionBegin;
12644482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
12654482741eSBarry Smith   PetscValidPointer(ksp,2);
1266988402f6SHong Zhang   if (((PetscObject)ts)->type_name == PETSC_NULL)
1267988402f6SHong Zhang     SETERRQ(PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
126829bbc08cSBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
126994b7f48cSBarry Smith   *ksp = ts->ksp;
1270d763cef2SBarry Smith   PetscFunctionReturn(0);
1271d763cef2SBarry Smith }
1272d763cef2SBarry Smith 
1273d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1274d763cef2SBarry Smith 
12754a2ae208SSatish Balay #undef __FUNCT__
1276adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1277adb62b0dSMatthew Knepley /*@
1278adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1279adb62b0dSMatthew Knepley    maximum time for iteration.
1280adb62b0dSMatthew Knepley 
1281adb62b0dSMatthew Knepley    Collective on TS
1282adb62b0dSMatthew Knepley 
1283adb62b0dSMatthew Knepley    Input Parameters:
1284adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1285adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1286adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1287adb62b0dSMatthew Knepley 
1288adb62b0dSMatthew Knepley    Level: intermediate
1289adb62b0dSMatthew Knepley 
1290adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1291adb62b0dSMatthew Knepley @*/
129263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1293adb62b0dSMatthew Knepley {
1294adb62b0dSMatthew Knepley   PetscFunctionBegin;
12954482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1296abc0a331SBarry Smith   if (maxsteps) {
12974482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1298adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1299adb62b0dSMatthew Knepley   }
1300abc0a331SBarry Smith   if (maxtime ) {
13014482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1302adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1303adb62b0dSMatthew Knepley   }
1304adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1305adb62b0dSMatthew Knepley }
1306adb62b0dSMatthew Knepley 
1307adb62b0dSMatthew Knepley #undef __FUNCT__
13084a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1309d763cef2SBarry Smith /*@
1310d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1311d763cef2SBarry Smith    maximum time for iteration.
1312d763cef2SBarry Smith 
1313d763cef2SBarry Smith    Collective on TS
1314d763cef2SBarry Smith 
1315d763cef2SBarry Smith    Input Parameters:
1316d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1317d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1318d763cef2SBarry Smith -  maxtime - final time to iterate to
1319d763cef2SBarry Smith 
1320d763cef2SBarry Smith    Options Database Keys:
1321d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1322d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1323d763cef2SBarry Smith 
1324d763cef2SBarry Smith    Notes:
1325d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1326d763cef2SBarry Smith 
1327d763cef2SBarry Smith    Level: intermediate
1328d763cef2SBarry Smith 
1329d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1330d763cef2SBarry Smith @*/
133163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1332d763cef2SBarry Smith {
1333d763cef2SBarry Smith   PetscFunctionBegin;
13344482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1335d763cef2SBarry Smith   ts->max_steps = maxsteps;
1336d763cef2SBarry Smith   ts->max_time  = maxtime;
1337d763cef2SBarry Smith   PetscFunctionReturn(0);
1338d763cef2SBarry Smith }
1339d763cef2SBarry Smith 
13404a2ae208SSatish Balay #undef __FUNCT__
13414a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1342d763cef2SBarry Smith /*@
1343d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1344d763cef2SBarry Smith    for use by the TS routines.
1345d763cef2SBarry Smith 
1346d763cef2SBarry Smith    Collective on TS and Vec
1347d763cef2SBarry Smith 
1348d763cef2SBarry Smith    Input Parameters:
1349d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1350d763cef2SBarry Smith -  x - the solution vector
1351d763cef2SBarry Smith 
1352d763cef2SBarry Smith    Level: beginner
1353d763cef2SBarry Smith 
1354d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1355d763cef2SBarry Smith @*/
135663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetSolution(TS ts,Vec x)
1357d763cef2SBarry Smith {
1358d763cef2SBarry Smith   PetscFunctionBegin;
13594482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
13604482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1361d763cef2SBarry Smith   ts->vec_sol        = ts->vec_sol_always = x;
1362d763cef2SBarry Smith   PetscFunctionReturn(0);
1363d763cef2SBarry Smith }
1364d763cef2SBarry Smith 
1365e74ef692SMatthew Knepley #undef __FUNCT__
1366e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1367ac226902SBarry Smith /*@C
1368000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
1369000e7ae3SMatthew Knepley   called once at the beginning of time stepping.
1370000e7ae3SMatthew Knepley 
1371000e7ae3SMatthew Knepley   Collective on TS
1372000e7ae3SMatthew Knepley 
1373000e7ae3SMatthew Knepley   Input Parameters:
1374000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1375000e7ae3SMatthew Knepley - func - The function
1376000e7ae3SMatthew Knepley 
1377000e7ae3SMatthew Knepley   Calling sequence of func:
1378000e7ae3SMatthew Knepley . func (TS ts);
1379000e7ae3SMatthew Knepley 
1380000e7ae3SMatthew Knepley   Level: intermediate
1381000e7ae3SMatthew Knepley 
1382000e7ae3SMatthew Knepley .keywords: TS, timestep
1383000e7ae3SMatthew Knepley @*/
138463dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1385000e7ae3SMatthew Knepley {
1386000e7ae3SMatthew Knepley   PetscFunctionBegin;
13874482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1388000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1389000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1390000e7ae3SMatthew Knepley }
1391000e7ae3SMatthew Knepley 
1392e74ef692SMatthew Knepley #undef __FUNCT__
1393e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPreStep"
1394000e7ae3SMatthew Knepley /*@
1395000e7ae3SMatthew Knepley   TSDefaultPreStep - The default pre-stepping function which does nothing.
1396000e7ae3SMatthew Knepley 
1397000e7ae3SMatthew Knepley   Collective on TS
1398000e7ae3SMatthew Knepley 
1399000e7ae3SMatthew Knepley   Input Parameters:
1400000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1401000e7ae3SMatthew Knepley 
1402000e7ae3SMatthew Knepley   Level: developer
1403000e7ae3SMatthew Knepley 
1404000e7ae3SMatthew Knepley .keywords: TS, timestep
1405000e7ae3SMatthew Knepley @*/
140663dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPreStep(TS ts)
1407000e7ae3SMatthew Knepley {
1408000e7ae3SMatthew Knepley   PetscFunctionBegin;
1409000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1410000e7ae3SMatthew Knepley }
1411000e7ae3SMatthew Knepley 
1412e74ef692SMatthew Knepley #undef __FUNCT__
1413e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1414ac226902SBarry Smith /*@C
1415000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
1416000e7ae3SMatthew Knepley   called once at the end of time stepping.
1417000e7ae3SMatthew Knepley 
1418000e7ae3SMatthew Knepley   Collective on TS
1419000e7ae3SMatthew Knepley 
1420000e7ae3SMatthew Knepley   Input Parameters:
1421000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1422000e7ae3SMatthew Knepley - func - The function
1423000e7ae3SMatthew Knepley 
1424000e7ae3SMatthew Knepley   Calling sequence of func:
1425000e7ae3SMatthew Knepley . func (TS ts);
1426000e7ae3SMatthew Knepley 
1427000e7ae3SMatthew Knepley   Level: intermediate
1428000e7ae3SMatthew Knepley 
1429000e7ae3SMatthew Knepley .keywords: TS, timestep
1430000e7ae3SMatthew Knepley @*/
143163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1432000e7ae3SMatthew Knepley {
1433000e7ae3SMatthew Knepley   PetscFunctionBegin;
14344482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1435000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1436000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1437000e7ae3SMatthew Knepley }
1438000e7ae3SMatthew Knepley 
1439e74ef692SMatthew Knepley #undef __FUNCT__
1440e74ef692SMatthew Knepley #define __FUNCT__ "TSDefaultPostStep"
1441000e7ae3SMatthew Knepley /*@
1442000e7ae3SMatthew Knepley   TSDefaultPostStep - The default post-stepping function which does nothing.
1443000e7ae3SMatthew Knepley 
1444000e7ae3SMatthew Knepley   Collective on TS
1445000e7ae3SMatthew Knepley 
1446000e7ae3SMatthew Knepley   Input Parameters:
1447000e7ae3SMatthew Knepley . ts  - The TS context obtained from TSCreate()
1448000e7ae3SMatthew Knepley 
1449000e7ae3SMatthew Knepley   Level: developer
1450000e7ae3SMatthew Knepley 
1451000e7ae3SMatthew Knepley .keywords: TS, timestep
1452000e7ae3SMatthew Knepley @*/
145363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSDefaultPostStep(TS ts)
1454000e7ae3SMatthew Knepley {
1455000e7ae3SMatthew Knepley   PetscFunctionBegin;
1456000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1457000e7ae3SMatthew Knepley }
1458000e7ae3SMatthew Knepley 
1459d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1460d763cef2SBarry Smith 
14614a2ae208SSatish Balay #undef __FUNCT__
1462a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet"
1463d763cef2SBarry Smith /*@C
1464a6570f20SBarry Smith    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1465d763cef2SBarry Smith    timestep to display the iteration's  progress.
1466d763cef2SBarry Smith 
1467d763cef2SBarry Smith    Collective on TS
1468d763cef2SBarry Smith 
1469d763cef2SBarry Smith    Input Parameters:
1470d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1471d763cef2SBarry Smith .  func - monitoring routine
1472329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1473b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1474b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1475b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1476d763cef2SBarry Smith 
1477d763cef2SBarry Smith    Calling sequence of func:
1478a7cc72afSBarry Smith $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1479d763cef2SBarry Smith 
1480d763cef2SBarry Smith +    ts - the TS context
1481d763cef2SBarry Smith .    steps - iteration number
14821f06c33eSBarry Smith .    time - current time
1483d763cef2SBarry Smith .    x - current iterate
1484d763cef2SBarry Smith -    mctx - [optional] monitoring context
1485d763cef2SBarry Smith 
1486d763cef2SBarry Smith    Notes:
1487d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1488d763cef2SBarry Smith    already has been loaded.
1489d763cef2SBarry Smith 
1490025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each TS object
1491025f1a04SBarry Smith 
1492d763cef2SBarry Smith    Level: intermediate
1493d763cef2SBarry Smith 
1494d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1495d763cef2SBarry Smith 
1496a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel()
1497d763cef2SBarry Smith @*/
1498a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void*))
1499d763cef2SBarry Smith {
1500d763cef2SBarry Smith   PetscFunctionBegin;
15014482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1502d763cef2SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) {
150329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1504d763cef2SBarry Smith   }
1505d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1506329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1507d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1508d763cef2SBarry Smith   PetscFunctionReturn(0);
1509d763cef2SBarry Smith }
1510d763cef2SBarry Smith 
15114a2ae208SSatish Balay #undef __FUNCT__
1512a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel"
1513d763cef2SBarry Smith /*@C
1514a6570f20SBarry Smith    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1515d763cef2SBarry Smith 
1516d763cef2SBarry Smith    Collective on TS
1517d763cef2SBarry Smith 
1518d763cef2SBarry Smith    Input Parameters:
1519d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1520d763cef2SBarry Smith 
1521d763cef2SBarry Smith    Notes:
1522d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1523d763cef2SBarry Smith 
1524d763cef2SBarry Smith    Level: intermediate
1525d763cef2SBarry Smith 
1526d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1527d763cef2SBarry Smith 
1528a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet()
1529d763cef2SBarry Smith @*/
1530a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorCancel(TS ts)
1531d763cef2SBarry Smith {
1532d952e501SBarry Smith   PetscErrorCode ierr;
1533d952e501SBarry Smith   PetscInt       i;
1534d952e501SBarry Smith 
1535d763cef2SBarry Smith   PetscFunctionBegin;
15364482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1537d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1538d952e501SBarry Smith     if (ts->mdestroy[i]) {
1539d952e501SBarry Smith       ierr = (*ts->mdestroy[i])(ts->monitorcontext[i]);CHKERRQ(ierr);
1540d952e501SBarry Smith     }
1541d952e501SBarry Smith   }
1542d763cef2SBarry Smith   ts->numbermonitors = 0;
1543d763cef2SBarry Smith   PetscFunctionReturn(0);
1544d763cef2SBarry Smith }
1545d763cef2SBarry Smith 
15464a2ae208SSatish Balay #undef __FUNCT__
1547a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault"
1548d8e5e3e6SSatish Balay /*@
1549a6570f20SBarry Smith    TSMonitorDefault - Sets the Default monitor
15505516499fSSatish Balay 
15515516499fSSatish Balay    Level: intermediate
155241251cbbSSatish Balay 
15535516499fSSatish Balay .keywords: TS, set, monitor
15545516499fSSatish Balay 
155541251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet()
155641251cbbSSatish Balay @*/
1557a6570f20SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
1558d763cef2SBarry Smith {
1559dfbe8321SBarry Smith   PetscErrorCode          ierr;
1560a34d58ebSBarry Smith   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor)ctx;
1561d132466eSBarry Smith 
1562d763cef2SBarry Smith   PetscFunctionBegin;
1563a34d58ebSBarry Smith   if (!ctx) {
15647adad957SLisandro Dalcin     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)ts)->comm,"stdout",0,&viewer);CHKERRQ(ierr);
1565a34d58ebSBarry Smith   }
1566a34d58ebSBarry Smith   ierr = PetscViewerASCIIMonitorPrintf(viewer,"timestep %D dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1567a34d58ebSBarry Smith   if (!ctx) {
1568a34d58ebSBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr);
1569a34d58ebSBarry Smith   }
1570d763cef2SBarry Smith   PetscFunctionReturn(0);
1571d763cef2SBarry Smith }
1572d763cef2SBarry Smith 
15734a2ae208SSatish Balay #undef __FUNCT__
15744a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1575d763cef2SBarry Smith /*@
1576d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1577d763cef2SBarry Smith 
1578d763cef2SBarry Smith    Collective on TS
1579d763cef2SBarry Smith 
1580d763cef2SBarry Smith    Input Parameter:
1581d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1582d763cef2SBarry Smith 
1583d763cef2SBarry Smith    Output Parameters:
1584d763cef2SBarry Smith +  steps - number of iterations until termination
1585142b95e3SSatish Balay -  ptime - time until termination
1586d763cef2SBarry Smith 
1587d763cef2SBarry Smith    Level: beginner
1588d763cef2SBarry Smith 
1589d763cef2SBarry Smith .keywords: TS, timestep, solve
1590d763cef2SBarry Smith 
1591d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1592d763cef2SBarry Smith @*/
159363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSStep(TS ts,PetscInt *steps,PetscReal *ptime)
1594d763cef2SBarry Smith {
1595dfbe8321SBarry Smith   PetscErrorCode ierr;
1596d763cef2SBarry Smith 
1597d763cef2SBarry Smith   PetscFunctionBegin;
15984482741eSBarry Smith   PetscValidHeaderSpecific(ts, TS_COOKIE,1);
1599d405a339SMatthew Knepley   if (!ts->setupcalled) {
1600d405a339SMatthew Knepley     ierr = TSSetUp(ts);CHKERRQ(ierr);
1601d405a339SMatthew Knepley   }
1602d405a339SMatthew Knepley 
1603d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1604d405a339SMatthew Knepley   ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
1605000e7ae3SMatthew Knepley   ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr);
1606d405a339SMatthew Knepley   ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
1607d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1608d405a339SMatthew Knepley 
16094bb05414SBarry Smith   if (!PetscPreLoadingOn) {
16107adad957SLisandro Dalcin     ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr);
1611d405a339SMatthew Knepley   }
1612d763cef2SBarry Smith   PetscFunctionReturn(0);
1613d763cef2SBarry Smith }
1614d763cef2SBarry Smith 
16154a2ae208SSatish Balay #undef __FUNCT__
16166a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve"
16176a4d4014SLisandro Dalcin /*@
16186a4d4014SLisandro Dalcin    TSSolve - Steps the requested number of timesteps.
16196a4d4014SLisandro Dalcin 
16206a4d4014SLisandro Dalcin    Collective on TS
16216a4d4014SLisandro Dalcin 
16226a4d4014SLisandro Dalcin    Input Parameter:
16236a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
16246a4d4014SLisandro Dalcin -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
16256a4d4014SLisandro Dalcin 
16266a4d4014SLisandro Dalcin    Level: beginner
16276a4d4014SLisandro Dalcin 
16286a4d4014SLisandro Dalcin .keywords: TS, timestep, solve
16296a4d4014SLisandro Dalcin 
16306a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep()
16316a4d4014SLisandro Dalcin @*/
16326a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSolve(TS ts, Vec x)
16336a4d4014SLisandro Dalcin {
16346a4d4014SLisandro Dalcin   PetscInt       steps;
16356a4d4014SLisandro Dalcin   PetscReal      ptime;
16366a4d4014SLisandro Dalcin   PetscErrorCode ierr;
16376a4d4014SLisandro Dalcin   PetscFunctionBegin;
16386a4d4014SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
16396a4d4014SLisandro Dalcin   /* set solution vector if provided */
16406a4d4014SLisandro Dalcin   if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); }
16416a4d4014SLisandro Dalcin   /* reset time step and iteration counters */
16426a4d4014SLisandro Dalcin   ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0;
16436a4d4014SLisandro Dalcin   /* steps the requested number of timesteps. */
16446a4d4014SLisandro Dalcin   ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr);
16456a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
16466a4d4014SLisandro Dalcin }
16476a4d4014SLisandro Dalcin 
16486a4d4014SLisandro Dalcin #undef __FUNCT__
16494a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1650d763cef2SBarry Smith /*
1651d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1652d763cef2SBarry Smith */
1653a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1654d763cef2SBarry Smith {
16556849ba73SBarry Smith   PetscErrorCode ierr;
1656a7cc72afSBarry Smith   PetscInt i,n = ts->numbermonitors;
1657d763cef2SBarry Smith 
1658d763cef2SBarry Smith   PetscFunctionBegin;
1659d763cef2SBarry Smith   for (i=0; i<n; i++) {
1660142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1661d763cef2SBarry Smith   }
1662d763cef2SBarry Smith   PetscFunctionReturn(0);
1663d763cef2SBarry Smith }
1664d763cef2SBarry Smith 
1665d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1666d763cef2SBarry Smith 
16674a2ae208SSatish Balay #undef __FUNCT__
1668a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate"
1669d763cef2SBarry Smith /*@C
1670a6570f20SBarry Smith    TSMonitorLGCreate - Creates a line graph context for use with
1671d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1672d763cef2SBarry Smith 
1673d763cef2SBarry Smith    Collective on TS
1674d763cef2SBarry Smith 
1675d763cef2SBarry Smith    Input Parameters:
1676d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1677d763cef2SBarry Smith .  label - the title to put in the title bar
16787c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1679d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1680d763cef2SBarry Smith 
1681d763cef2SBarry Smith    Output Parameter:
1682d763cef2SBarry Smith .  draw - the drawing context
1683d763cef2SBarry Smith 
1684d763cef2SBarry Smith    Options Database Key:
1685a6570f20SBarry Smith .  -ts_monitor_draw - automatically sets line graph monitor
1686d763cef2SBarry Smith 
1687d763cef2SBarry Smith    Notes:
1688a6570f20SBarry Smith    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1689d763cef2SBarry Smith 
1690d763cef2SBarry Smith    Level: intermediate
1691d763cef2SBarry Smith 
16927c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1693d763cef2SBarry Smith 
1694a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet()
16957c922b88SBarry Smith 
1696d763cef2SBarry Smith @*/
1697a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1698d763cef2SBarry Smith {
1699b0a32e0cSBarry Smith   PetscDraw      win;
1700dfbe8321SBarry Smith   PetscErrorCode ierr;
1701d763cef2SBarry Smith 
1702d763cef2SBarry Smith   PetscFunctionBegin;
1703b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1704b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1705b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1706b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1707d763cef2SBarry Smith 
170852e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1709d763cef2SBarry Smith   PetscFunctionReturn(0);
1710d763cef2SBarry Smith }
1711d763cef2SBarry Smith 
17124a2ae208SSatish Balay #undef __FUNCT__
1713a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG"
1714a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1715d763cef2SBarry Smith {
1716b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
171787828ca2SBarry Smith   PetscReal      x,y = ptime;
1718dfbe8321SBarry Smith   PetscErrorCode ierr;
1719d763cef2SBarry Smith 
1720d763cef2SBarry Smith   PetscFunctionBegin;
17217c922b88SBarry Smith   if (!monctx) {
17227c922b88SBarry Smith     MPI_Comm    comm;
1723b0a32e0cSBarry Smith     PetscViewer viewer;
17247c922b88SBarry Smith 
17257c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1726b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1727b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
17287c922b88SBarry Smith   }
17297c922b88SBarry Smith 
1730b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
173187828ca2SBarry Smith   x = (PetscReal)n;
1732b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1733d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1734b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1735d763cef2SBarry Smith   }
1736d763cef2SBarry Smith   PetscFunctionReturn(0);
1737d763cef2SBarry Smith }
1738d763cef2SBarry Smith 
17394a2ae208SSatish Balay #undef __FUNCT__
1740a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy"
1741d763cef2SBarry Smith /*@C
1742a6570f20SBarry Smith    TSMonitorLGDestroy - Destroys a line graph context that was created
1743a6570f20SBarry Smith    with TSMonitorLGCreate().
1744d763cef2SBarry Smith 
1745b0a32e0cSBarry Smith    Collective on PetscDrawLG
1746d763cef2SBarry Smith 
1747d763cef2SBarry Smith    Input Parameter:
1748d763cef2SBarry Smith .  draw - the drawing context
1749d763cef2SBarry Smith 
1750d763cef2SBarry Smith    Level: intermediate
1751d763cef2SBarry Smith 
1752d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1753d763cef2SBarry Smith 
1754a6570f20SBarry Smith .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1755d763cef2SBarry Smith @*/
1756a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorLGDestroy(PetscDrawLG drawlg)
1757d763cef2SBarry Smith {
1758b0a32e0cSBarry Smith   PetscDraw      draw;
1759dfbe8321SBarry Smith   PetscErrorCode ierr;
1760d763cef2SBarry Smith 
1761d763cef2SBarry Smith   PetscFunctionBegin;
1762b0a32e0cSBarry Smith   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
1763b0a32e0cSBarry Smith   ierr = PetscDrawDestroy(draw);CHKERRQ(ierr);
1764b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1765d763cef2SBarry Smith   PetscFunctionReturn(0);
1766d763cef2SBarry Smith }
1767d763cef2SBarry Smith 
17684a2ae208SSatish Balay #undef __FUNCT__
17694a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1770d763cef2SBarry Smith /*@
1771d763cef2SBarry Smith    TSGetTime - Gets the current time.
1772d763cef2SBarry Smith 
1773d763cef2SBarry Smith    Not Collective
1774d763cef2SBarry Smith 
1775d763cef2SBarry Smith    Input Parameter:
1776d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1777d763cef2SBarry Smith 
1778d763cef2SBarry Smith    Output Parameter:
1779d763cef2SBarry Smith .  t  - the current time
1780d763cef2SBarry Smith 
1781d763cef2SBarry Smith    Level: beginner
1782d763cef2SBarry Smith 
1783d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1784d763cef2SBarry Smith 
1785d763cef2SBarry Smith .keywords: TS, get, time
1786d763cef2SBarry Smith @*/
178763dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetTime(TS ts,PetscReal* t)
1788d763cef2SBarry Smith {
1789d763cef2SBarry Smith   PetscFunctionBegin;
17904482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
17914482741eSBarry Smith   PetscValidDoublePointer(t,2);
1792d763cef2SBarry Smith   *t = ts->ptime;
1793d763cef2SBarry Smith   PetscFunctionReturn(0);
1794d763cef2SBarry Smith }
1795d763cef2SBarry Smith 
17964a2ae208SSatish Balay #undef __FUNCT__
17976a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime"
17986a4d4014SLisandro Dalcin /*@
17996a4d4014SLisandro Dalcin    TSSetTime - Allows one to reset the time.
18006a4d4014SLisandro Dalcin 
18016a4d4014SLisandro Dalcin    Collective on TS
18026a4d4014SLisandro Dalcin 
18036a4d4014SLisandro Dalcin    Input Parameters:
18046a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
18056a4d4014SLisandro Dalcin -  time - the time
18066a4d4014SLisandro Dalcin 
18076a4d4014SLisandro Dalcin    Level: intermediate
18086a4d4014SLisandro Dalcin 
18096a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration()
18106a4d4014SLisandro Dalcin 
18116a4d4014SLisandro Dalcin .keywords: TS, set, time
18126a4d4014SLisandro Dalcin @*/
18136a4d4014SLisandro Dalcin PetscErrorCode PETSCTS_DLLEXPORT TSSetTime(TS ts, PetscReal t)
18146a4d4014SLisandro Dalcin {
18156a4d4014SLisandro Dalcin   PetscFunctionBegin;
18166a4d4014SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
18176a4d4014SLisandro Dalcin   ts->ptime = t;
18186a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
18196a4d4014SLisandro Dalcin }
18206a4d4014SLisandro Dalcin 
18216a4d4014SLisandro Dalcin #undef __FUNCT__
18224a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1823d763cef2SBarry Smith /*@C
1824d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1825d763cef2SBarry Smith    TS options in the database.
1826d763cef2SBarry Smith 
1827d763cef2SBarry Smith    Collective on TS
1828d763cef2SBarry Smith 
1829d763cef2SBarry Smith    Input Parameter:
1830d763cef2SBarry Smith +  ts     - The TS context
1831d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1832d763cef2SBarry Smith 
1833d763cef2SBarry Smith    Notes:
1834d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1835d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1836d763cef2SBarry Smith    hyphen.
1837d763cef2SBarry Smith 
1838d763cef2SBarry Smith    Level: advanced
1839d763cef2SBarry Smith 
1840d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1841d763cef2SBarry Smith 
1842d763cef2SBarry Smith .seealso: TSSetFromOptions()
1843d763cef2SBarry Smith 
1844d763cef2SBarry Smith @*/
184563dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSSetOptionsPrefix(TS ts,const char prefix[])
1846d763cef2SBarry Smith {
1847dfbe8321SBarry Smith   PetscErrorCode ierr;
1848d763cef2SBarry Smith 
1849d763cef2SBarry Smith   PetscFunctionBegin;
18504482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1851d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1852d763cef2SBarry Smith   switch(ts->problem_type) {
1853d763cef2SBarry Smith     case TS_NONLINEAR:
185459580b9cSBarry Smith       if (ts->snes) {
1855d763cef2SBarry Smith         ierr = SNESSetOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
185659580b9cSBarry Smith       }
1857d763cef2SBarry Smith       break;
1858d763cef2SBarry Smith     case TS_LINEAR:
185959580b9cSBarry Smith       if (ts->ksp) {
186094b7f48cSBarry Smith         ierr = KSPSetOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
186159580b9cSBarry Smith       }
1862d763cef2SBarry Smith       break;
1863d763cef2SBarry Smith   }
1864d763cef2SBarry Smith   PetscFunctionReturn(0);
1865d763cef2SBarry Smith }
1866d763cef2SBarry Smith 
1867d763cef2SBarry Smith 
18684a2ae208SSatish Balay #undef __FUNCT__
18694a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1870d763cef2SBarry Smith /*@C
1871d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1872d763cef2SBarry Smith    TS options in the database.
1873d763cef2SBarry Smith 
1874d763cef2SBarry Smith    Collective on TS
1875d763cef2SBarry Smith 
1876d763cef2SBarry Smith    Input Parameter:
1877d763cef2SBarry Smith +  ts     - The TS context
1878d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1879d763cef2SBarry Smith 
1880d763cef2SBarry Smith    Notes:
1881d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1882d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1883d763cef2SBarry Smith    hyphen.
1884d763cef2SBarry Smith 
1885d763cef2SBarry Smith    Level: advanced
1886d763cef2SBarry Smith 
1887d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1888d763cef2SBarry Smith 
1889d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1890d763cef2SBarry Smith 
1891d763cef2SBarry Smith @*/
189263dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSAppendOptionsPrefix(TS ts,const char prefix[])
1893d763cef2SBarry Smith {
1894dfbe8321SBarry Smith   PetscErrorCode ierr;
1895d763cef2SBarry Smith 
1896d763cef2SBarry Smith   PetscFunctionBegin;
18974482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
1898d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1899d763cef2SBarry Smith   switch(ts->problem_type) {
1900d763cef2SBarry Smith     case TS_NONLINEAR:
19011ac94b3bSBarry Smith       if (ts->snes) {
1902d763cef2SBarry Smith         ierr = SNESAppendOptionsPrefix(ts->snes,prefix);CHKERRQ(ierr);
19031ac94b3bSBarry Smith       }
1904d763cef2SBarry Smith       break;
1905d763cef2SBarry Smith     case TS_LINEAR:
19061ac94b3bSBarry Smith       if (ts->ksp) {
190794b7f48cSBarry Smith         ierr = KSPAppendOptionsPrefix(ts->ksp,prefix);CHKERRQ(ierr);
19081ac94b3bSBarry Smith       }
1909d763cef2SBarry Smith       break;
1910d763cef2SBarry Smith   }
1911d763cef2SBarry Smith   PetscFunctionReturn(0);
1912d763cef2SBarry Smith }
1913d763cef2SBarry Smith 
19144a2ae208SSatish Balay #undef __FUNCT__
19154a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1916d763cef2SBarry Smith /*@C
1917d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1918d763cef2SBarry Smith    TS options in the database.
1919d763cef2SBarry Smith 
1920d763cef2SBarry Smith    Not Collective
1921d763cef2SBarry Smith 
1922d763cef2SBarry Smith    Input Parameter:
1923d763cef2SBarry Smith .  ts - The TS context
1924d763cef2SBarry Smith 
1925d763cef2SBarry Smith    Output Parameter:
1926d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1927d763cef2SBarry Smith 
1928d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1929d763cef2SBarry Smith    sufficient length to hold the prefix.
1930d763cef2SBarry Smith 
1931d763cef2SBarry Smith    Level: intermediate
1932d763cef2SBarry Smith 
1933d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
1934d763cef2SBarry Smith 
1935d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
1936d763cef2SBarry Smith @*/
1937e060cb09SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSGetOptionsPrefix(TS ts,const char *prefix[])
1938d763cef2SBarry Smith {
1939dfbe8321SBarry Smith   PetscErrorCode ierr;
1940d763cef2SBarry Smith 
1941d763cef2SBarry Smith   PetscFunctionBegin;
19424482741eSBarry Smith   PetscValidHeaderSpecific(ts,TS_COOKIE,1);
19434482741eSBarry Smith   PetscValidPointer(prefix,2);
1944d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1945d763cef2SBarry Smith   PetscFunctionReturn(0);
1946d763cef2SBarry Smith }
1947d763cef2SBarry Smith 
19484a2ae208SSatish Balay #undef __FUNCT__
19494a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
1950d763cef2SBarry Smith /*@C
1951d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
1952d763cef2SBarry Smith 
1953d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
1954d763cef2SBarry Smith 
1955d763cef2SBarry Smith    Input Parameter:
1956d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
1957d763cef2SBarry Smith 
1958d763cef2SBarry Smith    Output Parameters:
1959d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
1960d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
1961d763cef2SBarry Smith - ctx - User-defined context for Jacobian evaluation routine
1962d763cef2SBarry Smith 
1963d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
1964d763cef2SBarry Smith 
1965d763cef2SBarry Smith    Level: intermediate
1966d763cef2SBarry Smith 
196726d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
1968d763cef2SBarry Smith 
1969d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
1970d763cef2SBarry Smith @*/
197163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSGetRHSJacobian(TS ts,Mat *J,Mat *M,void **ctx)
1972d763cef2SBarry Smith {
1973d763cef2SBarry Smith   PetscFunctionBegin;
197426d46c62SHong Zhang   if (J) *J = ts->Arhs;
197526d46c62SHong Zhang   if (M) *M = ts->B;
197626d46c62SHong Zhang   if (ctx) *ctx = ts->jacP;
1977d763cef2SBarry Smith   PetscFunctionReturn(0);
1978d763cef2SBarry Smith }
1979d763cef2SBarry Smith 
19801713a123SBarry Smith #undef __FUNCT__
1981a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution"
19821713a123SBarry Smith /*@C
1983a6570f20SBarry Smith    TSMonitorSolution - Monitors progress of the TS solvers by calling
19841713a123SBarry Smith    VecView() for the solution at each timestep
19851713a123SBarry Smith 
19861713a123SBarry Smith    Collective on TS
19871713a123SBarry Smith 
19881713a123SBarry Smith    Input Parameters:
19891713a123SBarry Smith +  ts - the TS context
19901713a123SBarry Smith .  step - current time-step
1991142b95e3SSatish Balay .  ptime - current time
19921713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
19931713a123SBarry Smith 
19941713a123SBarry Smith    Level: intermediate
19951713a123SBarry Smith 
19961713a123SBarry Smith .keywords: TS,  vector, monitor, view
19971713a123SBarry Smith 
1998a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
19991713a123SBarry Smith @*/
2000a6570f20SBarry Smith PetscErrorCode PETSCTS_DLLEXPORT TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
20011713a123SBarry Smith {
2002dfbe8321SBarry Smith   PetscErrorCode ierr;
20031713a123SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
20041713a123SBarry Smith 
20051713a123SBarry Smith   PetscFunctionBegin;
2006a34d58ebSBarry Smith   if (!dummy) {
20077adad957SLisandro Dalcin     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
20081713a123SBarry Smith   }
20091713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
20101713a123SBarry Smith   PetscFunctionReturn(0);
20111713a123SBarry Smith }
20121713a123SBarry Smith 
20131713a123SBarry Smith 
20141713a123SBarry Smith 
2015