xref: /petsc/src/ts/interface/ts.c (revision 3daf2de87434e83fb31fdc5ebc875393c291e8f9)
163dd3a1aSKris Buschelman 
2c6db04a5SJed Brown #include <private/tsimpl.h>        /*I "petscts.h"  I*/
3d763cef2SBarry Smith 
4d5ba7fb7SMatthew Knepley /* Logging support */
57087cfbeSBarry Smith PetscClassId  TS_CLASSID;
6166c7f25SBarry Smith PetscLogEvent  TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
7d405a339SMatthew Knepley 
84a2ae208SSatish Balay #undef __FUNCT__
9bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
10bdad233fSMatthew Knepley /*
11bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
12bdad233fSMatthew Knepley 
13bdad233fSMatthew Knepley   Collective on TS
14bdad233fSMatthew Knepley 
15bdad233fSMatthew Knepley   Input Parameter:
16bdad233fSMatthew Knepley . ts - The ts
17bdad233fSMatthew Knepley 
18bdad233fSMatthew Knepley   Level: intermediate
19bdad233fSMatthew Knepley 
20bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
21bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
22bdad233fSMatthew Knepley */
236849ba73SBarry Smith static PetscErrorCode TSSetTypeFromOptions(TS ts)
24bdad233fSMatthew Knepley {
25ace3abfcSBarry Smith   PetscBool      opt;
262fc52814SBarry Smith   const char     *defaultType;
27bdad233fSMatthew Knepley   char           typeName[256];
28dfbe8321SBarry Smith   PetscErrorCode ierr;
29bdad233fSMatthew Knepley 
30bdad233fSMatthew Knepley   PetscFunctionBegin;
317adad957SLisandro Dalcin   if (((PetscObject)ts)->type_name) {
327adad957SLisandro Dalcin     defaultType = ((PetscObject)ts)->type_name;
33bdad233fSMatthew Knepley   } else {
349596e0b4SJed Brown     defaultType = TSEULER;
35bdad233fSMatthew Knepley   }
36bdad233fSMatthew Knepley 
37cce0b1b2SLisandro Dalcin   if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
38bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
39a7cc72afSBarry Smith   if (opt) {
40bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
41bdad233fSMatthew Knepley   } else {
42bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
43bdad233fSMatthew Knepley   }
44bdad233fSMatthew Knepley   PetscFunctionReturn(0);
45bdad233fSMatthew Knepley }
46bdad233fSMatthew Knepley 
47bdad233fSMatthew Knepley #undef __FUNCT__
48bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
49bdad233fSMatthew Knepley /*@
50bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
51bdad233fSMatthew Knepley 
52bdad233fSMatthew Knepley    Collective on TS
53bdad233fSMatthew Knepley 
54bdad233fSMatthew Knepley    Input Parameter:
55bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
56bdad233fSMatthew Knepley 
57bdad233fSMatthew Knepley    Options Database Keys:
584d91e141SJed Brown +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP
59bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
60bdad233fSMatthew Knepley .  -ts_max_time time - maximum time to compute to
61bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
62bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
63a6570f20SBarry Smith -  -ts_monitor_draw - plot information at each timestep
64bdad233fSMatthew Knepley 
65bdad233fSMatthew Knepley    Level: beginner
66bdad233fSMatthew Knepley 
67bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
68bdad233fSMatthew Knepley 
69a313700dSBarry Smith .seealso: TSGetType()
70bdad233fSMatthew Knepley @*/
717087cfbeSBarry Smith PetscErrorCode  TSSetFromOptions(TS ts)
72bdad233fSMatthew Knepley {
73bdad233fSMatthew Knepley   PetscReal      dt;
74ace3abfcSBarry Smith   PetscBool      opt,flg;
75dfbe8321SBarry Smith   PetscErrorCode ierr;
76649052a6SBarry Smith   PetscViewer    monviewer;
77eabae89aSBarry Smith   char           monfilename[PETSC_MAX_PATH_LEN];
78089b2837SJed Brown   SNES           snes;
79bdad233fSMatthew Knepley 
80bdad233fSMatthew Knepley   PetscFunctionBegin;
810700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
827adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr);
83bdad233fSMatthew Knepley 
84bdad233fSMatthew Knepley     /* Handle generic TS options */
85bdad233fSMatthew Knepley     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
86bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
87bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr);
88bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr);
89*3daf2de8SJed Brown     if (opt) {ierr = TSSetInitialTimeStep(ts,ts->ptime,dt);CHKERRQ(ierr);}
90193ac0bcSJed Brown     ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr);
91193ac0bcSJed Brown     ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections","",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr);
92193ac0bcSJed Brown     ierr = PetscOptionsBool("-ts_error_if_step_failed","Error if no step succeeds","",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr);
93bdad233fSMatthew Knepley 
94bdad233fSMatthew Knepley     /* Monitor options */
95a6570f20SBarry Smith     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
96eabae89aSBarry Smith     if (flg) {
97649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr);
98649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
99bdad233fSMatthew Knepley     }
10090d69ab7SBarry Smith     opt  = PETSC_FALSE;
101acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
102a7cc72afSBarry Smith     if (opt) {
103a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
104bdad233fSMatthew Knepley     }
10590d69ab7SBarry Smith     opt  = PETSC_FALSE;
106acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
107a7cc72afSBarry Smith     if (opt) {
108a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
109bdad233fSMatthew Knepley     }
110bdad233fSMatthew Knepley 
111bdad233fSMatthew Knepley     /* Handle TS type options */
112bdad233fSMatthew Knepley     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
113bdad233fSMatthew Knepley 
114bdad233fSMatthew Knepley     /* Handle specific TS options */
115abc0a331SBarry Smith     if (ts->ops->setfromoptions) {
116bdad233fSMatthew Knepley       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
117bdad233fSMatthew Knepley     }
1185d973c19SBarry Smith 
1195d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
1205d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
121bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
122bdad233fSMatthew Knepley 
123089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
124bdad233fSMatthew Knepley   /* Handle subobject options */
125089b2837SJed Brown   if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
126089b2837SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
127bdad233fSMatthew Knepley   PetscFunctionReturn(0);
128bdad233fSMatthew Knepley }
129bdad233fSMatthew Knepley 
130bdad233fSMatthew Knepley #undef  __FUNCT__
131bdad233fSMatthew Knepley #define __FUNCT__ "TSViewFromOptions"
132bdad233fSMatthew Knepley /*@
133bdad233fSMatthew Knepley   TSViewFromOptions - This function visualizes the ts based upon user options.
134bdad233fSMatthew Knepley 
135bdad233fSMatthew Knepley   Collective on TS
136bdad233fSMatthew Knepley 
137bdad233fSMatthew Knepley   Input Parameter:
138bdad233fSMatthew Knepley . ts - The ts
139bdad233fSMatthew Knepley 
140bdad233fSMatthew Knepley   Level: intermediate
141bdad233fSMatthew Knepley 
142bdad233fSMatthew Knepley .keywords: TS, view, options, database
143bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSView()
144bdad233fSMatthew Knepley @*/
1457087cfbeSBarry Smith PetscErrorCode  TSViewFromOptions(TS ts,const char title[])
146bdad233fSMatthew Knepley {
147bdad233fSMatthew Knepley   PetscViewer    viewer;
148bdad233fSMatthew Knepley   PetscDraw      draw;
149ace3abfcSBarry Smith   PetscBool      opt = PETSC_FALSE;
150e10c49a3SBarry Smith   char           fileName[PETSC_MAX_PATH_LEN];
151dfbe8321SBarry Smith   PetscErrorCode ierr;
152bdad233fSMatthew Knepley 
153bdad233fSMatthew Knepley   PetscFunctionBegin;
1547adad957SLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr);
155eabae89aSBarry Smith   if (opt && !PetscPreLoadingOn) {
1567adad957SLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr);
157bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
1586bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
159bdad233fSMatthew Knepley   }
1608e83347fSKai Germaschewski   opt = PETSC_FALSE;
161acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr);
162a7cc72afSBarry Smith   if (opt) {
1637adad957SLisandro Dalcin     ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr);
164bdad233fSMatthew Knepley     ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr);
165a7cc72afSBarry Smith     if (title) {
1661836bdbcSSatish Balay       ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr);
167bdad233fSMatthew Knepley     } else {
168bdad233fSMatthew Knepley       ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr);
1697adad957SLisandro Dalcin       ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr);
170bdad233fSMatthew Knepley     }
171bdad233fSMatthew Knepley     ierr = TSView(ts, viewer);CHKERRQ(ierr);
172bdad233fSMatthew Knepley     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
173bdad233fSMatthew Knepley     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
1746bf464f9SBarry Smith     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
175bdad233fSMatthew Knepley   }
176bdad233fSMatthew Knepley   PetscFunctionReturn(0);
177bdad233fSMatthew Knepley }
178bdad233fSMatthew Knepley 
179bdad233fSMatthew Knepley #undef __FUNCT__
1804a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
181a7a1495cSBarry Smith /*@
1828c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
183a7a1495cSBarry Smith       set with TSSetRHSJacobian().
184a7a1495cSBarry Smith 
185a7a1495cSBarry Smith    Collective on TS and Vec
186a7a1495cSBarry Smith 
187a7a1495cSBarry Smith    Input Parameters:
188316643e7SJed Brown +  ts - the TS context
189a7a1495cSBarry Smith .  t - current timestep
190a7a1495cSBarry Smith -  x - input vector
191a7a1495cSBarry Smith 
192a7a1495cSBarry Smith    Output Parameters:
193a7a1495cSBarry Smith +  A - Jacobian matrix
194a7a1495cSBarry Smith .  B - optional preconditioning matrix
195a7a1495cSBarry Smith -  flag - flag indicating matrix structure
196a7a1495cSBarry Smith 
197a7a1495cSBarry Smith    Notes:
198a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
199a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
200a7a1495cSBarry Smith 
20194b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
202a7a1495cSBarry Smith    flag parameter.
203a7a1495cSBarry Smith 
204a7a1495cSBarry Smith    Level: developer
205a7a1495cSBarry Smith 
206a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
207a7a1495cSBarry Smith 
20894b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
209a7a1495cSBarry Smith @*/
2107087cfbeSBarry Smith PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
211a7a1495cSBarry Smith {
212dfbe8321SBarry Smith   PetscErrorCode ierr;
213a7a1495cSBarry Smith 
214a7a1495cSBarry Smith   PetscFunctionBegin;
2150700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2160700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
217c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
218000e7ae3SMatthew Knepley   if (ts->ops->rhsjacobian) {
219d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
220a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
221a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
222000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
223a7a1495cSBarry Smith     PetscStackPop;
224d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
225a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2260700a824SBarry Smith     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
2270700a824SBarry Smith     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
228ef66eb69SBarry Smith   } else {
229214bc6a2SJed Brown     ierr = MatZeroEntries(*A);CHKERRQ(ierr);
230214bc6a2SJed Brown     if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);}
231214bc6a2SJed Brown     *flg = SAME_NONZERO_PATTERN;
232ef66eb69SBarry Smith   }
233a7a1495cSBarry Smith   PetscFunctionReturn(0);
234a7a1495cSBarry Smith }
235a7a1495cSBarry Smith 
2364a2ae208SSatish Balay #undef __FUNCT__
2374a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
238316643e7SJed Brown /*@
239d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
240d763cef2SBarry Smith 
241316643e7SJed Brown    Collective on TS and Vec
242316643e7SJed Brown 
243316643e7SJed Brown    Input Parameters:
244316643e7SJed Brown +  ts - the TS context
245316643e7SJed Brown .  t - current time
246316643e7SJed Brown -  x - state vector
247316643e7SJed Brown 
248316643e7SJed Brown    Output Parameter:
249316643e7SJed Brown .  y - right hand side
250316643e7SJed Brown 
251316643e7SJed Brown    Note:
252316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
253316643e7SJed Brown    is used internally within the nonlinear solvers.
254316643e7SJed Brown 
255316643e7SJed Brown    Level: developer
256316643e7SJed Brown 
257316643e7SJed Brown .keywords: TS, compute
258316643e7SJed Brown 
259316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction()
260316643e7SJed Brown @*/
261dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
262d763cef2SBarry Smith {
263dfbe8321SBarry Smith   PetscErrorCode ierr;
264d763cef2SBarry Smith 
265d763cef2SBarry Smith   PetscFunctionBegin;
2660700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2670700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2680700a824SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
269d763cef2SBarry Smith 
270d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
271000e7ae3SMatthew Knepley   if (ts->ops->rhsfunction) {
272d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
273000e7ae3SMatthew Knepley     ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
274d763cef2SBarry Smith     PetscStackPop;
275214bc6a2SJed Brown   } else {
276214bc6a2SJed Brown     ierr = VecZeroEntries(y);CHKERRQ(ierr);
277b2cd27e8SJed Brown   }
278d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
279d763cef2SBarry Smith   PetscFunctionReturn(0);
280d763cef2SBarry Smith }
281d763cef2SBarry Smith 
2824a2ae208SSatish Balay #undef __FUNCT__
283214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSVec_Private"
284214bc6a2SJed Brown static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
285214bc6a2SJed Brown {
286214bc6a2SJed Brown   PetscErrorCode ierr;
287214bc6a2SJed Brown 
288214bc6a2SJed Brown   PetscFunctionBegin;
289214bc6a2SJed Brown   if (!ts->Frhs) {
290214bc6a2SJed Brown     ierr = VecDuplicate(ts->vec_sol,&ts->Frhs);CHKERRQ(ierr);
291214bc6a2SJed Brown   }
292214bc6a2SJed Brown   *Frhs = ts->Frhs;
293214bc6a2SJed Brown   PetscFunctionReturn(0);
294214bc6a2SJed Brown }
295214bc6a2SJed Brown 
296214bc6a2SJed Brown #undef __FUNCT__
297214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSMats_Private"
298214bc6a2SJed Brown static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
299214bc6a2SJed Brown {
300214bc6a2SJed Brown   PetscErrorCode ierr;
301214bc6a2SJed Brown   Mat A,B;
302214bc6a2SJed Brown 
303214bc6a2SJed Brown   PetscFunctionBegin;
304214bc6a2SJed Brown   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
305214bc6a2SJed Brown   if (Arhs) {
306214bc6a2SJed Brown     if (!ts->Arhs) {
307214bc6a2SJed Brown       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
308214bc6a2SJed Brown     }
309214bc6a2SJed Brown     *Arhs = ts->Arhs;
310214bc6a2SJed Brown   }
311214bc6a2SJed Brown   if (Brhs) {
312214bc6a2SJed Brown     if (!ts->Brhs) {
313214bc6a2SJed Brown       ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
314214bc6a2SJed Brown     }
315214bc6a2SJed Brown     *Brhs = ts->Brhs;
316214bc6a2SJed Brown   }
317214bc6a2SJed Brown   PetscFunctionReturn(0);
318214bc6a2SJed Brown }
319214bc6a2SJed Brown 
320214bc6a2SJed Brown #undef __FUNCT__
321316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction"
322316643e7SJed Brown /*@
323316643e7SJed Brown    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
324316643e7SJed Brown 
325316643e7SJed Brown    Collective on TS and Vec
326316643e7SJed Brown 
327316643e7SJed Brown    Input Parameters:
328316643e7SJed Brown +  ts - the TS context
329316643e7SJed Brown .  t - current time
330316643e7SJed Brown .  X - state vector
331214bc6a2SJed Brown .  Xdot - time derivative of state vector
332214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
333316643e7SJed Brown 
334316643e7SJed Brown    Output Parameter:
335316643e7SJed Brown .  Y - right hand side
336316643e7SJed Brown 
337316643e7SJed Brown    Note:
338316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
339316643e7SJed Brown    is used internally within the nonlinear solvers.
340316643e7SJed Brown 
341316643e7SJed Brown    If the user did did not write their equations in implicit form, this
342316643e7SJed Brown    function recasts them in implicit form.
343316643e7SJed Brown 
344316643e7SJed Brown    Level: developer
345316643e7SJed Brown 
346316643e7SJed Brown .keywords: TS, compute
347316643e7SJed Brown 
348316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction()
349316643e7SJed Brown @*/
350214bc6a2SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
351316643e7SJed Brown {
352316643e7SJed Brown   PetscErrorCode ierr;
353316643e7SJed Brown 
354316643e7SJed Brown   PetscFunctionBegin;
3550700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3560700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
3570700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
3580700a824SBarry Smith   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
359316643e7SJed Brown 
360316643e7SJed Brown   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
361316643e7SJed Brown   if (ts->ops->ifunction) {
362316643e7SJed Brown     PetscStackPush("TS user implicit function");
363316643e7SJed Brown     ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr);
364316643e7SJed Brown     PetscStackPop;
365214bc6a2SJed Brown   }
366214bc6a2SJed Brown   if (imex) {
367214bc6a2SJed Brown     if (!ts->ops->ifunction) {ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);}
368316643e7SJed Brown   } else {
369214bc6a2SJed Brown     if (!ts->ops->ifunction) {
370089b2837SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr);
371ace68cafSJed Brown       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
372214bc6a2SJed Brown     } else {
373214bc6a2SJed Brown       Vec Frhs;
374214bc6a2SJed Brown       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
375214bc6a2SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr);
376214bc6a2SJed Brown       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
377316643e7SJed Brown     }
3784a6899ffSJed Brown   }
379316643e7SJed Brown   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
380316643e7SJed Brown   PetscFunctionReturn(0);
381316643e7SJed Brown }
382316643e7SJed Brown 
383316643e7SJed Brown #undef __FUNCT__
384316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian"
385316643e7SJed Brown /*@
386316643e7SJed Brown    TSComputeIJacobian - Evaluates the Jacobian of the DAE
387316643e7SJed Brown 
388316643e7SJed Brown    Collective on TS and Vec
389316643e7SJed Brown 
390316643e7SJed Brown    Input
391316643e7SJed Brown       Input Parameters:
392316643e7SJed Brown +  ts - the TS context
393316643e7SJed Brown .  t - current timestep
394316643e7SJed Brown .  X - state vector
395316643e7SJed Brown .  Xdot - time derivative of state vector
396214bc6a2SJed Brown .  shift - shift to apply, see note below
397214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
398316643e7SJed Brown 
399316643e7SJed Brown    Output Parameters:
400316643e7SJed Brown +  A - Jacobian matrix
401316643e7SJed Brown .  B - optional preconditioning matrix
402316643e7SJed Brown -  flag - flag indicating matrix structure
403316643e7SJed Brown 
404316643e7SJed Brown    Notes:
405316643e7SJed Brown    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
406316643e7SJed Brown 
407316643e7SJed Brown    dF/dX + shift*dF/dXdot
408316643e7SJed Brown 
409316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
410316643e7SJed Brown    is used internally within the nonlinear solvers.
411316643e7SJed Brown 
412316643e7SJed Brown    Level: developer
413316643e7SJed Brown 
414316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix
415316643e7SJed Brown 
416316643e7SJed Brown .seealso:  TSSetIJacobian()
41763495f91SJed Brown @*/
418214bc6a2SJed Brown PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
419316643e7SJed Brown {
420316643e7SJed Brown   PetscErrorCode ierr;
421316643e7SJed Brown 
422316643e7SJed Brown   PetscFunctionBegin;
4230700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4240700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
4250700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
426316643e7SJed Brown   PetscValidPointer(A,6);
4270700a824SBarry Smith   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
428316643e7SJed Brown   PetscValidPointer(B,7);
4290700a824SBarry Smith   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
430316643e7SJed Brown   PetscValidPointer(flg,8);
431316643e7SJed Brown 
4324e684422SJed Brown   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
433316643e7SJed Brown   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
434316643e7SJed Brown   if (ts->ops->ijacobian) {
435316643e7SJed Brown     PetscStackPush("TS user implicit Jacobian");
436316643e7SJed Brown     ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr);
437316643e7SJed Brown     PetscStackPop;
438214bc6a2SJed Brown     /* make sure user returned a correct Jacobian and preconditioner */
439214bc6a2SJed Brown     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
440214bc6a2SJed Brown     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
4414a6899ffSJed Brown   }
442214bc6a2SJed Brown   if (imex) {
443214bc6a2SJed Brown     if (!ts->ops->ijacobian) {  /* system was written as Xdot = F(t,X) */
444214bc6a2SJed Brown       ierr = MatZeroEntries(*A);CHKERRQ(ierr);
445214bc6a2SJed Brown       ierr = MatShift(*A,1.0);CHKERRQ(ierr);
446214bc6a2SJed Brown       if (*A != *B) {
447214bc6a2SJed Brown         ierr = MatZeroEntries(*B);CHKERRQ(ierr);
448214bc6a2SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
449214bc6a2SJed Brown       }
450214bc6a2SJed Brown       *flg = SAME_PRECONDITIONER;
451214bc6a2SJed Brown     }
452214bc6a2SJed Brown   } else {
453214bc6a2SJed Brown     if (!ts->ops->ijacobian) {
454214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr);
455214bc6a2SJed Brown       ierr = MatScale(*A,-1);CHKERRQ(ierr);
456214bc6a2SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
457316643e7SJed Brown       if (*A != *B) {
458316643e7SJed Brown         ierr = MatScale(*B,-1);CHKERRQ(ierr);
459316643e7SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
460316643e7SJed Brown       }
461214bc6a2SJed Brown     } else if (ts->ops->rhsjacobian) {
462214bc6a2SJed Brown       Mat Arhs,Brhs;
463214bc6a2SJed Brown       MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
464214bc6a2SJed Brown       ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
465214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
466214bc6a2SJed Brown       axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
467214bc6a2SJed Brown       ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr);
468214bc6a2SJed Brown       if (*A != *B) {
469214bc6a2SJed Brown         ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr);
470214bc6a2SJed Brown       }
471214bc6a2SJed Brown       *flg = PetscMin(*flg,flg2);
472214bc6a2SJed Brown     }
473316643e7SJed Brown   }
474316643e7SJed Brown   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
475316643e7SJed Brown   PetscFunctionReturn(0);
476316643e7SJed Brown }
477316643e7SJed Brown 
478316643e7SJed Brown #undef __FUNCT__
4794a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
480d763cef2SBarry Smith /*@C
481d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
482d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
483d763cef2SBarry Smith 
4843f9fe445SBarry Smith     Logically Collective on TS
485d763cef2SBarry Smith 
486d763cef2SBarry Smith     Input Parameters:
487d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
488ca94891dSJed Brown .   r - vector to put the computed right hand side (or PETSC_NULL to have it created)
489d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
490d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
491d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
492d763cef2SBarry Smith 
493d763cef2SBarry Smith     Calling sequence of func:
49487828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
495d763cef2SBarry Smith 
496d763cef2SBarry Smith +   t - current timestep
497d763cef2SBarry Smith .   u - input vector
498d763cef2SBarry Smith .   F - function vector
499d763cef2SBarry Smith -   ctx - [optional] user-defined function context
500d763cef2SBarry Smith 
501d763cef2SBarry Smith     Important:
50276f2fa84SHong Zhang     The user MUST call either this routine or TSSetMatrices().
503d763cef2SBarry Smith 
504d763cef2SBarry Smith     Level: beginner
505d763cef2SBarry Smith 
506d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
507d763cef2SBarry Smith 
50876f2fa84SHong Zhang .seealso: TSSetMatrices()
509d763cef2SBarry Smith @*/
510089b2837SJed Brown PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
511d763cef2SBarry Smith {
512089b2837SJed Brown   PetscErrorCode ierr;
513089b2837SJed Brown   SNES           snes;
514d763cef2SBarry Smith 
515089b2837SJed Brown   PetscFunctionBegin;
5160700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
517ca94891dSJed Brown   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
5184e684422SJed Brown   if (f)   ts->ops->rhsfunction = f;
5194e684422SJed Brown   if (ctx) ts->funP             = ctx;
520089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
521089b2837SJed Brown   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
522d763cef2SBarry Smith   PetscFunctionReturn(0);
523d763cef2SBarry Smith }
524d763cef2SBarry Smith 
5254a2ae208SSatish Balay #undef __FUNCT__
5264a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
527d763cef2SBarry Smith /*@C
528d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
529d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
53076f2fa84SHong Zhang    Use TSSetMatrices() for linear problems.
531d763cef2SBarry Smith 
5323f9fe445SBarry Smith    Logically Collective on TS
533d763cef2SBarry Smith 
534d763cef2SBarry Smith    Input Parameters:
535d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
536d763cef2SBarry Smith .  A   - Jacobian matrix
537d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
538d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
539d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
540d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
541d763cef2SBarry Smith 
542d763cef2SBarry Smith    Calling sequence of func:
54387828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
544d763cef2SBarry Smith 
545d763cef2SBarry Smith +  t - current timestep
546d763cef2SBarry Smith .  u - input vector
547d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
548d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
549d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
55094b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
551d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
552d763cef2SBarry Smith 
553d763cef2SBarry Smith    Notes:
55494b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
555d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
556d763cef2SBarry Smith 
557d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
558d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
55956335db2SHong Zhang    completely new matrix structure (not just different matrix elements)
560d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
561d763cef2SBarry Smith    throughout the global iterations.
562d763cef2SBarry Smith 
563d763cef2SBarry Smith    Level: beginner
564d763cef2SBarry Smith 
565d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
566d763cef2SBarry Smith 
567d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
56876f2fa84SHong Zhang           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices()
569d763cef2SBarry Smith 
570d763cef2SBarry Smith @*/
571089b2837SJed Brown PetscErrorCode  TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
572d763cef2SBarry Smith {
573277b19d0SLisandro Dalcin   PetscErrorCode ierr;
574089b2837SJed Brown   SNES           snes;
575277b19d0SLisandro Dalcin 
576d763cef2SBarry Smith   PetscFunctionBegin;
5770700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
578277b19d0SLisandro Dalcin   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
579277b19d0SLisandro Dalcin   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
580277b19d0SLisandro Dalcin   if (A) PetscCheckSameComm(ts,1,A,2);
581277b19d0SLisandro Dalcin   if (B) PetscCheckSameComm(ts,1,B,3);
582d763cef2SBarry Smith 
583277b19d0SLisandro Dalcin   if (f)   ts->ops->rhsjacobian = f;
584277b19d0SLisandro Dalcin   if (ctx) ts->jacP             = ctx;
585089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
586089b2837SJed Brown   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
587d763cef2SBarry Smith   PetscFunctionReturn(0);
588d763cef2SBarry Smith }
589d763cef2SBarry Smith 
590316643e7SJed Brown 
591316643e7SJed Brown #undef __FUNCT__
592316643e7SJed Brown #define __FUNCT__ "TSSetIFunction"
593316643e7SJed Brown /*@C
594316643e7SJed Brown    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
595316643e7SJed Brown 
5963f9fe445SBarry Smith    Logically Collective on TS
597316643e7SJed Brown 
598316643e7SJed Brown    Input Parameters:
599316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
600ca94891dSJed Brown .  r   - vector to hold the residual (or PETSC_NULL to have it created internally)
601316643e7SJed Brown .  f   - the function evaluation routine
602316643e7SJed Brown -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
603316643e7SJed Brown 
604316643e7SJed Brown    Calling sequence of f:
605316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
606316643e7SJed Brown 
607316643e7SJed Brown +  t   - time at step/stage being solved
608316643e7SJed Brown .  u   - state vector
609316643e7SJed Brown .  u_t - time derivative of state vector
610316643e7SJed Brown .  F   - function vector
611316643e7SJed Brown -  ctx - [optional] user-defined context for matrix evaluation routine
612316643e7SJed Brown 
613316643e7SJed Brown    Important:
614316643e7SJed Brown    The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices().  This routine must be used when not solving an ODE.
615316643e7SJed Brown 
616316643e7SJed Brown    Level: beginner
617316643e7SJed Brown 
618316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian
619316643e7SJed Brown 
620316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian()
621316643e7SJed Brown @*/
622089b2837SJed Brown PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
623316643e7SJed Brown {
624089b2837SJed Brown   PetscErrorCode ierr;
625089b2837SJed Brown   SNES           snes;
626316643e7SJed Brown 
627316643e7SJed Brown   PetscFunctionBegin;
6280700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
629ca94891dSJed Brown   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
630089b2837SJed Brown   if (f)   ts->ops->ifunction = f;
631089b2837SJed Brown   if (ctx) ts->funP           = ctx;
632089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
633089b2837SJed Brown   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
634089b2837SJed Brown   PetscFunctionReturn(0);
635089b2837SJed Brown }
636089b2837SJed Brown 
637089b2837SJed Brown #undef __FUNCT__
638089b2837SJed Brown #define __FUNCT__ "TSGetIFunction"
639089b2837SJed Brown /*@C
640089b2837SJed Brown    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
641089b2837SJed Brown 
642089b2837SJed Brown    Not Collective
643089b2837SJed Brown 
644089b2837SJed Brown    Input Parameter:
645089b2837SJed Brown .  ts - the TS context
646089b2837SJed Brown 
647089b2837SJed Brown    Output Parameter:
648089b2837SJed Brown +  r - vector to hold residual (or PETSC_NULL)
649089b2837SJed Brown .  func - the function to compute residual (or PETSC_NULL)
650089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
651089b2837SJed Brown 
652089b2837SJed Brown    Level: advanced
653089b2837SJed Brown 
654089b2837SJed Brown .keywords: TS, nonlinear, get, function
655089b2837SJed Brown 
656089b2837SJed Brown .seealso: TSSetIFunction(), SNESGetFunction()
657089b2837SJed Brown @*/
658089b2837SJed Brown PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
659089b2837SJed Brown {
660089b2837SJed Brown   PetscErrorCode ierr;
661089b2837SJed Brown   SNES snes;
662089b2837SJed Brown 
663089b2837SJed Brown   PetscFunctionBegin;
664089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
665089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
666089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
667089b2837SJed Brown   if (func) *func = ts->ops->ifunction;
668089b2837SJed Brown   if (ctx)  *ctx  = ts->funP;
669089b2837SJed Brown   PetscFunctionReturn(0);
670089b2837SJed Brown }
671089b2837SJed Brown 
672089b2837SJed Brown #undef __FUNCT__
673089b2837SJed Brown #define __FUNCT__ "TSGetRHSFunction"
674089b2837SJed Brown /*@C
675089b2837SJed Brown    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
676089b2837SJed Brown 
677089b2837SJed Brown    Not Collective
678089b2837SJed Brown 
679089b2837SJed Brown    Input Parameter:
680089b2837SJed Brown .  ts - the TS context
681089b2837SJed Brown 
682089b2837SJed Brown    Output Parameter:
683089b2837SJed Brown +  r - vector to hold computed right hand side (or PETSC_NULL)
684089b2837SJed Brown .  func - the function to compute right hand side (or PETSC_NULL)
685089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
686089b2837SJed Brown 
687089b2837SJed Brown    Level: advanced
688089b2837SJed Brown 
689089b2837SJed Brown .keywords: TS, nonlinear, get, function
690089b2837SJed Brown 
691089b2837SJed Brown .seealso: TSSetRhsfunction(), SNESGetFunction()
692089b2837SJed Brown @*/
693089b2837SJed Brown PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
694089b2837SJed Brown {
695089b2837SJed Brown   PetscErrorCode ierr;
696089b2837SJed Brown   SNES snes;
697089b2837SJed Brown 
698089b2837SJed Brown   PetscFunctionBegin;
699089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
700089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
701089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
702089b2837SJed Brown   if (func) *func = ts->ops->rhsfunction;
703089b2837SJed Brown   if (ctx)  *ctx  = ts->funP;
704316643e7SJed Brown   PetscFunctionReturn(0);
705316643e7SJed Brown }
706316643e7SJed Brown 
707316643e7SJed Brown #undef __FUNCT__
708316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian"
709316643e7SJed Brown /*@C
710a4f0a591SBarry Smith    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
711a4f0a591SBarry Smith         you provided with TSSetIFunction().
712316643e7SJed Brown 
7133f9fe445SBarry Smith    Logically Collective on TS
714316643e7SJed Brown 
715316643e7SJed Brown    Input Parameters:
716316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
717316643e7SJed Brown .  A   - Jacobian matrix
718316643e7SJed Brown .  B   - preconditioning matrix for A (may be same as A)
719316643e7SJed Brown .  f   - the Jacobian evaluation routine
720316643e7SJed Brown -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
721316643e7SJed Brown 
722316643e7SJed Brown    Calling sequence of f:
7231b4a444bSJed Brown $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
724316643e7SJed Brown 
725316643e7SJed Brown +  t    - time at step/stage being solved
7261b4a444bSJed Brown .  U    - state vector
7271b4a444bSJed Brown .  U_t  - time derivative of state vector
728316643e7SJed Brown .  a    - shift
7291b4a444bSJed Brown .  A    - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
730316643e7SJed Brown .  B    - preconditioning matrix for A, may be same as A
731316643e7SJed Brown .  flag - flag indicating information about the preconditioner matrix
732316643e7SJed Brown           structure (same as flag in KSPSetOperators())
733316643e7SJed Brown -  ctx  - [optional] user-defined context for matrix evaluation routine
734316643e7SJed Brown 
735316643e7SJed Brown    Notes:
736316643e7SJed Brown    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
737316643e7SJed Brown 
738a4f0a591SBarry Smith    The matrix dF/dU + a*dF/dU_t you provide turns out to be
739a4f0a591SBarry Smith    the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved.
740a4f0a591SBarry Smith    The time integrator internally approximates U_t by W+a*U where the positive "shift"
741a4f0a591SBarry Smith    a and vector W depend on the integration method, step size, and past states. For example with
742a4f0a591SBarry Smith    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
743a4f0a591SBarry Smith    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
744a4f0a591SBarry Smith 
745316643e7SJed Brown    Level: beginner
746316643e7SJed Brown 
747316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian
748316643e7SJed Brown 
749316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian()
750316643e7SJed Brown 
751316643e7SJed Brown @*/
7527087cfbeSBarry Smith PetscErrorCode  TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
753316643e7SJed Brown {
754316643e7SJed Brown   PetscErrorCode ierr;
755089b2837SJed Brown   SNES           snes;
756316643e7SJed Brown 
757316643e7SJed Brown   PetscFunctionBegin;
7580700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7590700a824SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
7600700a824SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
761316643e7SJed Brown   if (A) PetscCheckSameComm(ts,1,A,2);
762316643e7SJed Brown   if (B) PetscCheckSameComm(ts,1,B,3);
763316643e7SJed Brown   if (f)   ts->ops->ijacobian = f;
764316643e7SJed Brown   if (ctx) ts->jacP           = ctx;
765089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
766089b2837SJed Brown   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
767316643e7SJed Brown   PetscFunctionReturn(0);
768316643e7SJed Brown }
769316643e7SJed Brown 
7704a2ae208SSatish Balay #undef __FUNCT__
7714a2ae208SSatish Balay #define __FUNCT__ "TSView"
7727e2c5f70SBarry Smith /*@C
773d763cef2SBarry Smith     TSView - Prints the TS data structure.
774d763cef2SBarry Smith 
7754c49b128SBarry Smith     Collective on TS
776d763cef2SBarry Smith 
777d763cef2SBarry Smith     Input Parameters:
778d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
779d763cef2SBarry Smith -   viewer - visualization context
780d763cef2SBarry Smith 
781d763cef2SBarry Smith     Options Database Key:
782d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
783d763cef2SBarry Smith 
784d763cef2SBarry Smith     Notes:
785d763cef2SBarry Smith     The available visualization contexts include
786b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
787b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
788d763cef2SBarry Smith          output where only the first processor opens
789d763cef2SBarry Smith          the file.  All other processors send their
790d763cef2SBarry Smith          data to the first processor to print.
791d763cef2SBarry Smith 
792d763cef2SBarry Smith     The user can open an alternative visualization context with
793b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
794d763cef2SBarry Smith 
795d763cef2SBarry Smith     Level: beginner
796d763cef2SBarry Smith 
797d763cef2SBarry Smith .keywords: TS, timestep, view
798d763cef2SBarry Smith 
799b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
800d763cef2SBarry Smith @*/
8017087cfbeSBarry Smith PetscErrorCode  TSView(TS ts,PetscViewer viewer)
802d763cef2SBarry Smith {
803dfbe8321SBarry Smith   PetscErrorCode ierr;
804a313700dSBarry Smith   const TSType   type;
805089b2837SJed Brown   PetscBool      iascii,isstring,isundials;
806d763cef2SBarry Smith 
807d763cef2SBarry Smith   PetscFunctionBegin;
8080700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8093050cee2SBarry Smith   if (!viewer) {
8107adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
8113050cee2SBarry Smith   }
8120700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
813c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
814fd16b177SBarry Smith 
8152692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8162692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
81732077d6dSBarry Smith   if (iascii) {
818317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr);
819000e7ae3SMatthew Knepley     if (ts->ops->view) {
820b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
821000e7ae3SMatthew Knepley       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
822b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
823d763cef2SBarry Smith     }
82477431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
825a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
826d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
82777431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
828193ac0bcSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->max_snes_failures);CHKERRQ(ierr);
829d763cef2SBarry Smith     }
83077431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
831193ac0bcSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
8320f5bd95cSBarry Smith   } else if (isstring) {
833a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
834b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
835d763cef2SBarry Smith   }
836b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
837089b2837SJed Brown   ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
838089b2837SJed Brown   if (!isundials && ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
839b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
840d763cef2SBarry Smith   PetscFunctionReturn(0);
841d763cef2SBarry Smith }
842d763cef2SBarry Smith 
843d763cef2SBarry Smith 
8444a2ae208SSatish Balay #undef __FUNCT__
8454a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
846b07ff414SBarry Smith /*@
847d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
848d763cef2SBarry Smith    the timesteppers.
849d763cef2SBarry Smith 
8503f9fe445SBarry Smith    Logically Collective on TS
851d763cef2SBarry Smith 
852d763cef2SBarry Smith    Input Parameters:
853d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
854d763cef2SBarry Smith -  usrP - optional user context
855d763cef2SBarry Smith 
856d763cef2SBarry Smith    Level: intermediate
857d763cef2SBarry Smith 
858d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
859d763cef2SBarry Smith 
860d763cef2SBarry Smith .seealso: TSGetApplicationContext()
861d763cef2SBarry Smith @*/
8627087cfbeSBarry Smith PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
863d763cef2SBarry Smith {
864d763cef2SBarry Smith   PetscFunctionBegin;
8650700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
866d763cef2SBarry Smith   ts->user = usrP;
867d763cef2SBarry Smith   PetscFunctionReturn(0);
868d763cef2SBarry Smith }
869d763cef2SBarry Smith 
8704a2ae208SSatish Balay #undef __FUNCT__
8714a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
872b07ff414SBarry Smith /*@
873d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
874d763cef2SBarry Smith     timestepper.
875d763cef2SBarry Smith 
876d763cef2SBarry Smith     Not Collective
877d763cef2SBarry Smith 
878d763cef2SBarry Smith     Input Parameter:
879d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
880d763cef2SBarry Smith 
881d763cef2SBarry Smith     Output Parameter:
882d763cef2SBarry Smith .   usrP - user context
883d763cef2SBarry Smith 
884d763cef2SBarry Smith     Level: intermediate
885d763cef2SBarry Smith 
886d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
887d763cef2SBarry Smith 
888d763cef2SBarry Smith .seealso: TSSetApplicationContext()
889d763cef2SBarry Smith @*/
890e71120c6SJed Brown PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
891d763cef2SBarry Smith {
892d763cef2SBarry Smith   PetscFunctionBegin;
8930700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
894e71120c6SJed Brown   *(void**)usrP = ts->user;
895d763cef2SBarry Smith   PetscFunctionReturn(0);
896d763cef2SBarry Smith }
897d763cef2SBarry Smith 
8984a2ae208SSatish Balay #undef __FUNCT__
8994a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
900d763cef2SBarry Smith /*@
901d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
902d763cef2SBarry Smith 
903d763cef2SBarry Smith    Not Collective
904d763cef2SBarry Smith 
905d763cef2SBarry Smith    Input Parameter:
906d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
907d763cef2SBarry Smith 
908d763cef2SBarry Smith    Output Parameter:
909d763cef2SBarry Smith .  iter - number steps so far
910d763cef2SBarry Smith 
911d763cef2SBarry Smith    Level: intermediate
912d763cef2SBarry Smith 
913d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
914d763cef2SBarry Smith @*/
9157087cfbeSBarry Smith PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
916d763cef2SBarry Smith {
917d763cef2SBarry Smith   PetscFunctionBegin;
9180700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
9194482741eSBarry Smith   PetscValidIntPointer(iter,2);
920d763cef2SBarry Smith   *iter = ts->steps;
921d763cef2SBarry Smith   PetscFunctionReturn(0);
922d763cef2SBarry Smith }
923d763cef2SBarry Smith 
9244a2ae208SSatish Balay #undef __FUNCT__
9254a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
926d763cef2SBarry Smith /*@
927d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
928d763cef2SBarry Smith    as well as the initial time.
929d763cef2SBarry Smith 
9303f9fe445SBarry Smith    Logically Collective on TS
931d763cef2SBarry Smith 
932d763cef2SBarry Smith    Input Parameters:
933d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
934d763cef2SBarry Smith .  initial_time - the initial time
935d763cef2SBarry Smith -  time_step - the size of the timestep
936d763cef2SBarry Smith 
937d763cef2SBarry Smith    Level: intermediate
938d763cef2SBarry Smith 
939d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
940d763cef2SBarry Smith 
941d763cef2SBarry Smith .keywords: TS, set, initial, timestep
942d763cef2SBarry Smith @*/
9437087cfbeSBarry Smith PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
944d763cef2SBarry Smith {
945e144a568SJed Brown   PetscErrorCode ierr;
946e144a568SJed Brown 
947d763cef2SBarry Smith   PetscFunctionBegin;
9480700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
949e144a568SJed Brown   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
950d763cef2SBarry Smith   ts->initial_time_step = time_step;
951d763cef2SBarry Smith   ts->ptime             = initial_time;
952d763cef2SBarry Smith   PetscFunctionReturn(0);
953d763cef2SBarry Smith }
954d763cef2SBarry Smith 
9554a2ae208SSatish Balay #undef __FUNCT__
9564a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
957d763cef2SBarry Smith /*@
958d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
959d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
960d763cef2SBarry Smith 
9613f9fe445SBarry Smith    Logically Collective on TS
962d763cef2SBarry Smith 
963d763cef2SBarry Smith    Input Parameters:
964d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
965d763cef2SBarry Smith -  time_step - the size of the timestep
966d763cef2SBarry Smith 
967d763cef2SBarry Smith    Level: intermediate
968d763cef2SBarry Smith 
969d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
970d763cef2SBarry Smith 
971d763cef2SBarry Smith .keywords: TS, set, timestep
972d763cef2SBarry Smith @*/
9737087cfbeSBarry Smith PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
974d763cef2SBarry Smith {
975d763cef2SBarry Smith   PetscFunctionBegin;
9760700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
977c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,time_step,2);
978d763cef2SBarry Smith   ts->time_step      = time_step;
979193ac0bcSJed Brown   ts->next_time_step = time_step;
980d763cef2SBarry Smith   PetscFunctionReturn(0);
981d763cef2SBarry Smith }
982d763cef2SBarry Smith 
9834a2ae208SSatish Balay #undef __FUNCT__
9844a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
985d763cef2SBarry Smith /*@
986d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
987d763cef2SBarry Smith 
988d763cef2SBarry Smith    Not Collective
989d763cef2SBarry Smith 
990d763cef2SBarry Smith    Input Parameter:
991d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
992d763cef2SBarry Smith 
993d763cef2SBarry Smith    Output Parameter:
994d763cef2SBarry Smith .  dt - the current timestep size
995d763cef2SBarry Smith 
996d763cef2SBarry Smith    Level: intermediate
997d763cef2SBarry Smith 
998d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
999d763cef2SBarry Smith 
1000d763cef2SBarry Smith .keywords: TS, get, timestep
1001d763cef2SBarry Smith @*/
10027087cfbeSBarry Smith PetscErrorCode  TSGetTimeStep(TS ts,PetscReal* dt)
1003d763cef2SBarry Smith {
1004d763cef2SBarry Smith   PetscFunctionBegin;
10050700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10064482741eSBarry Smith   PetscValidDoublePointer(dt,2);
1007d763cef2SBarry Smith   *dt = ts->time_step;
1008d763cef2SBarry Smith   PetscFunctionReturn(0);
1009d763cef2SBarry Smith }
1010d763cef2SBarry Smith 
10114a2ae208SSatish Balay #undef __FUNCT__
10124a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
1013d8e5e3e6SSatish Balay /*@
1014d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
1015d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
1016d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
1017d763cef2SBarry Smith    the solution at the next timestep has been calculated.
1018d763cef2SBarry Smith 
1019d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
1020d763cef2SBarry Smith 
1021d763cef2SBarry Smith    Input Parameter:
1022d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1023d763cef2SBarry Smith 
1024d763cef2SBarry Smith    Output Parameter:
1025d763cef2SBarry Smith .  v - the vector containing the solution
1026d763cef2SBarry Smith 
1027d763cef2SBarry Smith    Level: intermediate
1028d763cef2SBarry Smith 
1029d763cef2SBarry Smith .seealso: TSGetTimeStep()
1030d763cef2SBarry Smith 
1031d763cef2SBarry Smith .keywords: TS, timestep, get, solution
1032d763cef2SBarry Smith @*/
10337087cfbeSBarry Smith PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1034d763cef2SBarry Smith {
1035d763cef2SBarry Smith   PetscFunctionBegin;
10360700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10374482741eSBarry Smith   PetscValidPointer(v,2);
10388737fe31SLisandro Dalcin   *v = ts->vec_sol;
1039d763cef2SBarry Smith   PetscFunctionReturn(0);
1040d763cef2SBarry Smith }
1041d763cef2SBarry Smith 
1042bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
10434a2ae208SSatish Balay #undef __FUNCT__
1044bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
1045d8e5e3e6SSatish Balay /*@
1046bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
1047d763cef2SBarry Smith 
1048bdad233fSMatthew Knepley   Not collective
1049d763cef2SBarry Smith 
1050bdad233fSMatthew Knepley   Input Parameters:
1051bdad233fSMatthew Knepley + ts   - The TS
1052bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1053d763cef2SBarry Smith .vb
1054d763cef2SBarry Smith          U_t = A U
1055d763cef2SBarry Smith          U_t = A(t) U
1056d763cef2SBarry Smith          U_t = F(t,U)
1057d763cef2SBarry Smith .ve
1058d763cef2SBarry Smith 
1059d763cef2SBarry Smith    Level: beginner
1060d763cef2SBarry Smith 
1061bdad233fSMatthew Knepley .keywords: TS, problem type
1062bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1063d763cef2SBarry Smith @*/
10647087cfbeSBarry Smith PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1065a7cc72afSBarry Smith {
10669e2a6581SJed Brown   PetscErrorCode ierr;
10679e2a6581SJed Brown 
1068d763cef2SBarry Smith   PetscFunctionBegin;
10690700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1070bdad233fSMatthew Knepley   ts->problem_type = type;
10719e2a6581SJed Brown   if (type == TS_LINEAR) {
10729e2a6581SJed Brown     SNES snes;
10739e2a6581SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
10749e2a6581SJed Brown     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
10759e2a6581SJed Brown   }
1076d763cef2SBarry Smith   PetscFunctionReturn(0);
1077d763cef2SBarry Smith }
1078d763cef2SBarry Smith 
1079bdad233fSMatthew Knepley #undef __FUNCT__
1080bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
1081bdad233fSMatthew Knepley /*@C
1082bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
1083bdad233fSMatthew Knepley 
1084bdad233fSMatthew Knepley   Not collective
1085bdad233fSMatthew Knepley 
1086bdad233fSMatthew Knepley   Input Parameter:
1087bdad233fSMatthew Knepley . ts   - The TS
1088bdad233fSMatthew Knepley 
1089bdad233fSMatthew Knepley   Output Parameter:
1090bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1091bdad233fSMatthew Knepley .vb
1092089b2837SJed Brown          M U_t = A U
1093089b2837SJed Brown          M(t) U_t = A(t) U
1094bdad233fSMatthew Knepley          U_t = F(t,U)
1095bdad233fSMatthew Knepley .ve
1096bdad233fSMatthew Knepley 
1097bdad233fSMatthew Knepley    Level: beginner
1098bdad233fSMatthew Knepley 
1099bdad233fSMatthew Knepley .keywords: TS, problem type
1100bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1101bdad233fSMatthew Knepley @*/
11027087cfbeSBarry Smith PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1103a7cc72afSBarry Smith {
1104bdad233fSMatthew Knepley   PetscFunctionBegin;
11050700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
11064482741eSBarry Smith   PetscValidIntPointer(type,2);
1107bdad233fSMatthew Knepley   *type = ts->problem_type;
1108bdad233fSMatthew Knepley   PetscFunctionReturn(0);
1109bdad233fSMatthew Knepley }
1110d763cef2SBarry Smith 
11114a2ae208SSatish Balay #undef __FUNCT__
11124a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
1113d763cef2SBarry Smith /*@
1114d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
1115d763cef2SBarry Smith    of a timestepper.
1116d763cef2SBarry Smith 
1117d763cef2SBarry Smith    Collective on TS
1118d763cef2SBarry Smith 
1119d763cef2SBarry Smith    Input Parameter:
1120d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1121d763cef2SBarry Smith 
1122d763cef2SBarry Smith    Notes:
1123d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
1124d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
1125d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
1126d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
1127d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
1128d763cef2SBarry Smith 
1129d763cef2SBarry Smith    Level: advanced
1130d763cef2SBarry Smith 
1131d763cef2SBarry Smith .keywords: TS, timestep, setup
1132d763cef2SBarry Smith 
1133d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
1134d763cef2SBarry Smith @*/
11357087cfbeSBarry Smith PetscErrorCode  TSSetUp(TS ts)
1136d763cef2SBarry Smith {
1137dfbe8321SBarry Smith   PetscErrorCode ierr;
1138d763cef2SBarry Smith 
1139d763cef2SBarry Smith   PetscFunctionBegin;
11400700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1141277b19d0SLisandro Dalcin   if (ts->setupcalled) PetscFunctionReturn(0);
1142277b19d0SLisandro Dalcin 
11437adad957SLisandro Dalcin   if (!((PetscObject)ts)->type_name) {
11449596e0b4SJed Brown     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1145d763cef2SBarry Smith   }
1146277b19d0SLisandro Dalcin 
1147277b19d0SLisandro Dalcin   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1148277b19d0SLisandro Dalcin 
1149277b19d0SLisandro Dalcin   if (ts->ops->setup) {
1150000e7ae3SMatthew Knepley     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1151277b19d0SLisandro Dalcin   }
1152277b19d0SLisandro Dalcin 
1153277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_TRUE;
1154277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1155277b19d0SLisandro Dalcin }
1156277b19d0SLisandro Dalcin 
1157277b19d0SLisandro Dalcin #undef __FUNCT__
1158277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset"
1159277b19d0SLisandro Dalcin /*@
1160277b19d0SLisandro Dalcin    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1161277b19d0SLisandro Dalcin 
1162277b19d0SLisandro Dalcin    Collective on TS
1163277b19d0SLisandro Dalcin 
1164277b19d0SLisandro Dalcin    Input Parameter:
1165277b19d0SLisandro Dalcin .  ts - the TS context obtained from TSCreate()
1166277b19d0SLisandro Dalcin 
1167277b19d0SLisandro Dalcin    Level: beginner
1168277b19d0SLisandro Dalcin 
1169277b19d0SLisandro Dalcin .keywords: TS, timestep, reset
1170277b19d0SLisandro Dalcin 
1171277b19d0SLisandro Dalcin .seealso: TSCreate(), TSSetup(), TSDestroy()
1172277b19d0SLisandro Dalcin @*/
1173277b19d0SLisandro Dalcin PetscErrorCode  TSReset(TS ts)
1174277b19d0SLisandro Dalcin {
1175277b19d0SLisandro Dalcin   PetscErrorCode ierr;
1176277b19d0SLisandro Dalcin 
1177277b19d0SLisandro Dalcin   PetscFunctionBegin;
1178277b19d0SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1179277b19d0SLisandro Dalcin   if (ts->ops->reset) {
1180277b19d0SLisandro Dalcin     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1181277b19d0SLisandro Dalcin   }
1182277b19d0SLisandro Dalcin   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
11834e684422SJed Brown   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
11844e684422SJed Brown   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1185214bc6a2SJed Brown   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
11866bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1187277b19d0SLisandro Dalcin   if (ts->work) {ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);}
1188277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_FALSE;
1189d763cef2SBarry Smith   PetscFunctionReturn(0);
1190d763cef2SBarry Smith }
1191d763cef2SBarry Smith 
11924a2ae208SSatish Balay #undef __FUNCT__
11934a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
1194d8e5e3e6SSatish Balay /*@
1195d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
1196d763cef2SBarry Smith    with TSCreate().
1197d763cef2SBarry Smith 
1198d763cef2SBarry Smith    Collective on TS
1199d763cef2SBarry Smith 
1200d763cef2SBarry Smith    Input Parameter:
1201d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1202d763cef2SBarry Smith 
1203d763cef2SBarry Smith    Level: beginner
1204d763cef2SBarry Smith 
1205d763cef2SBarry Smith .keywords: TS, timestepper, destroy
1206d763cef2SBarry Smith 
1207d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
1208d763cef2SBarry Smith @*/
12096bf464f9SBarry Smith PetscErrorCode  TSDestroy(TS *ts)
1210d763cef2SBarry Smith {
12116849ba73SBarry Smith   PetscErrorCode ierr;
1212d763cef2SBarry Smith 
1213d763cef2SBarry Smith   PetscFunctionBegin;
12146bf464f9SBarry Smith   if (!*ts) PetscFunctionReturn(0);
12156bf464f9SBarry Smith   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
12166bf464f9SBarry Smith   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1217d763cef2SBarry Smith 
12186bf464f9SBarry Smith   ierr = TSReset((*ts));CHKERRQ(ierr);
1219277b19d0SLisandro Dalcin 
1220be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
12216bf464f9SBarry Smith   ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr);
12226bf464f9SBarry Smith   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
12236d4c513bSLisandro Dalcin 
12246bf464f9SBarry Smith   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
12256bf464f9SBarry Smith   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
12266bf464f9SBarry Smith   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
12276d4c513bSLisandro Dalcin 
1228a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1229d763cef2SBarry Smith   PetscFunctionReturn(0);
1230d763cef2SBarry Smith }
1231d763cef2SBarry Smith 
12324a2ae208SSatish Balay #undef __FUNCT__
12334a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
1234d8e5e3e6SSatish Balay /*@
1235d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1236d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
1237d763cef2SBarry Smith 
1238d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
1239d763cef2SBarry Smith 
1240d763cef2SBarry Smith    Input Parameter:
1241d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1242d763cef2SBarry Smith 
1243d763cef2SBarry Smith    Output Parameter:
1244d763cef2SBarry Smith .  snes - the nonlinear solver context
1245d763cef2SBarry Smith 
1246d763cef2SBarry Smith    Notes:
1247d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
1248d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
124994b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
1250d763cef2SBarry Smith 
1251d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
1252d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
1253d763cef2SBarry Smith 
1254d763cef2SBarry Smith    Level: beginner
1255d763cef2SBarry Smith 
1256d763cef2SBarry Smith .keywords: timestep, get, SNES
1257d763cef2SBarry Smith @*/
12587087cfbeSBarry Smith PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1259d763cef2SBarry Smith {
1260d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1261d372ba47SLisandro Dalcin 
1262d763cef2SBarry Smith   PetscFunctionBegin;
12630700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12644482741eSBarry Smith   PetscValidPointer(snes,2);
1265d372ba47SLisandro Dalcin   if (!ts->snes) {
1266d372ba47SLisandro Dalcin     ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
1267d372ba47SLisandro Dalcin     ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr);
1268d372ba47SLisandro Dalcin     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
12699e2a6581SJed Brown     if (ts->problem_type == TS_LINEAR) {
12709e2a6581SJed Brown       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
12719e2a6581SJed Brown     }
1272d372ba47SLisandro Dalcin   }
1273d763cef2SBarry Smith   *snes = ts->snes;
1274d763cef2SBarry Smith   PetscFunctionReturn(0);
1275d763cef2SBarry Smith }
1276d763cef2SBarry Smith 
12774a2ae208SSatish Balay #undef __FUNCT__
127894b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
1279d8e5e3e6SSatish Balay /*@
128094b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
1281d763cef2SBarry Smith    a TS (timestepper) context.
1282d763cef2SBarry Smith 
128394b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1284d763cef2SBarry Smith 
1285d763cef2SBarry Smith    Input Parameter:
1286d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1287d763cef2SBarry Smith 
1288d763cef2SBarry Smith    Output Parameter:
128994b7f48cSBarry Smith .  ksp - the nonlinear solver context
1290d763cef2SBarry Smith 
1291d763cef2SBarry Smith    Notes:
129294b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1293d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1294d763cef2SBarry Smith    KSP and PC contexts as well.
1295d763cef2SBarry Smith 
129694b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
129794b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1298d763cef2SBarry Smith 
1299d763cef2SBarry Smith    Level: beginner
1300d763cef2SBarry Smith 
130194b7f48cSBarry Smith .keywords: timestep, get, KSP
1302d763cef2SBarry Smith @*/
13037087cfbeSBarry Smith PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1304d763cef2SBarry Smith {
1305d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1306089b2837SJed Brown   SNES           snes;
1307d372ba47SLisandro Dalcin 
1308d763cef2SBarry Smith   PetscFunctionBegin;
13090700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13104482741eSBarry Smith   PetscValidPointer(ksp,2);
131117186662SBarry Smith   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1312e32f2f54SBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1313089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1314089b2837SJed Brown   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
1315d763cef2SBarry Smith   PetscFunctionReturn(0);
1316d763cef2SBarry Smith }
1317d763cef2SBarry Smith 
1318d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1319d763cef2SBarry Smith 
13204a2ae208SSatish Balay #undef __FUNCT__
1321adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1322adb62b0dSMatthew Knepley /*@
1323adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1324adb62b0dSMatthew Knepley    maximum time for iteration.
1325adb62b0dSMatthew Knepley 
13263f9fe445SBarry Smith    Not Collective
1327adb62b0dSMatthew Knepley 
1328adb62b0dSMatthew Knepley    Input Parameters:
1329adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1330adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1331adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1332adb62b0dSMatthew Knepley 
1333adb62b0dSMatthew Knepley    Level: intermediate
1334adb62b0dSMatthew Knepley 
1335adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1336adb62b0dSMatthew Knepley @*/
13377087cfbeSBarry Smith PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1338adb62b0dSMatthew Knepley {
1339adb62b0dSMatthew Knepley   PetscFunctionBegin;
13400700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1341abc0a331SBarry Smith   if (maxsteps) {
13424482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1343adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1344adb62b0dSMatthew Knepley   }
1345abc0a331SBarry Smith   if (maxtime ) {
13464482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1347adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1348adb62b0dSMatthew Knepley   }
1349adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1350adb62b0dSMatthew Knepley }
1351adb62b0dSMatthew Knepley 
1352adb62b0dSMatthew Knepley #undef __FUNCT__
13534a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1354d763cef2SBarry Smith /*@
1355d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1356d763cef2SBarry Smith    maximum time for iteration.
1357d763cef2SBarry Smith 
13583f9fe445SBarry Smith    Logically Collective on TS
1359d763cef2SBarry Smith 
1360d763cef2SBarry Smith    Input Parameters:
1361d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1362d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1363d763cef2SBarry Smith -  maxtime - final time to iterate to
1364d763cef2SBarry Smith 
1365d763cef2SBarry Smith    Options Database Keys:
1366d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1367d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1368d763cef2SBarry Smith 
1369d763cef2SBarry Smith    Notes:
1370d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1371d763cef2SBarry Smith 
1372d763cef2SBarry Smith    Level: intermediate
1373d763cef2SBarry Smith 
1374d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1375d763cef2SBarry Smith @*/
13767087cfbeSBarry Smith PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1377d763cef2SBarry Smith {
1378d763cef2SBarry Smith   PetscFunctionBegin;
13790700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1380c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1381c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,maxtime,2);
1382d763cef2SBarry Smith   ts->max_steps = maxsteps;
1383d763cef2SBarry Smith   ts->max_time  = maxtime;
1384d763cef2SBarry Smith   PetscFunctionReturn(0);
1385d763cef2SBarry Smith }
1386d763cef2SBarry Smith 
13874a2ae208SSatish Balay #undef __FUNCT__
13884a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1389d763cef2SBarry Smith /*@
1390d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1391d763cef2SBarry Smith    for use by the TS routines.
1392d763cef2SBarry Smith 
13933f9fe445SBarry Smith    Logically Collective on TS and Vec
1394d763cef2SBarry Smith 
1395d763cef2SBarry Smith    Input Parameters:
1396d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1397d763cef2SBarry Smith -  x - the solution vector
1398d763cef2SBarry Smith 
1399d763cef2SBarry Smith    Level: beginner
1400d763cef2SBarry Smith 
1401d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1402d763cef2SBarry Smith @*/
14037087cfbeSBarry Smith PetscErrorCode  TSSetSolution(TS ts,Vec x)
1404d763cef2SBarry Smith {
14058737fe31SLisandro Dalcin   PetscErrorCode ierr;
14068737fe31SLisandro Dalcin 
1407d763cef2SBarry Smith   PetscFunctionBegin;
14080700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
14090700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
14108737fe31SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
14116bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
14128737fe31SLisandro Dalcin   ts->vec_sol = x;
1413d763cef2SBarry Smith   PetscFunctionReturn(0);
1414d763cef2SBarry Smith }
1415d763cef2SBarry Smith 
1416e74ef692SMatthew Knepley #undef __FUNCT__
1417e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1418ac226902SBarry Smith /*@C
1419000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
14203f2090d5SJed Brown   called once at the beginning of each time step.
1421000e7ae3SMatthew Knepley 
14223f9fe445SBarry Smith   Logically Collective on TS
1423000e7ae3SMatthew Knepley 
1424000e7ae3SMatthew Knepley   Input Parameters:
1425000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1426000e7ae3SMatthew Knepley - func - The function
1427000e7ae3SMatthew Knepley 
1428000e7ae3SMatthew Knepley   Calling sequence of func:
1429000e7ae3SMatthew Knepley . func (TS ts);
1430000e7ae3SMatthew Knepley 
1431000e7ae3SMatthew Knepley   Level: intermediate
1432000e7ae3SMatthew Knepley 
1433000e7ae3SMatthew Knepley .keywords: TS, timestep
1434000e7ae3SMatthew Knepley @*/
14357087cfbeSBarry Smith PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1436000e7ae3SMatthew Knepley {
1437000e7ae3SMatthew Knepley   PetscFunctionBegin;
14380700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1439000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1440000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1441000e7ae3SMatthew Knepley }
1442000e7ae3SMatthew Knepley 
1443e74ef692SMatthew Knepley #undef __FUNCT__
14443f2090d5SJed Brown #define __FUNCT__ "TSPreStep"
14453f2090d5SJed Brown /*@C
14463f2090d5SJed Brown   TSPreStep - Runs the user-defined pre-step function.
14473f2090d5SJed Brown 
14483f2090d5SJed Brown   Collective on TS
14493f2090d5SJed Brown 
14503f2090d5SJed Brown   Input Parameters:
14513f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
14523f2090d5SJed Brown 
14533f2090d5SJed Brown   Notes:
14543f2090d5SJed Brown   TSPreStep() is typically used within time stepping implementations,
14553f2090d5SJed Brown   so most users would not generally call this routine themselves.
14563f2090d5SJed Brown 
14573f2090d5SJed Brown   Level: developer
14583f2090d5SJed Brown 
14593f2090d5SJed Brown .keywords: TS, timestep
14603f2090d5SJed Brown @*/
14617087cfbeSBarry Smith PetscErrorCode  TSPreStep(TS ts)
14623f2090d5SJed Brown {
14633f2090d5SJed Brown   PetscErrorCode ierr;
14643f2090d5SJed Brown 
14653f2090d5SJed Brown   PetscFunctionBegin;
14660700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
146772ac3e02SJed Brown   if (ts->ops->prestep) {
14683f2090d5SJed Brown     PetscStackPush("TS PreStep function");
14693f2090d5SJed Brown     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
14703f2090d5SJed Brown     PetscStackPop;
1471312ce896SJed Brown   }
14723f2090d5SJed Brown   PetscFunctionReturn(0);
14733f2090d5SJed Brown }
14743f2090d5SJed Brown 
14753f2090d5SJed Brown #undef __FUNCT__
1476e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1477ac226902SBarry Smith /*@C
1478000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
14793f2090d5SJed Brown   called once at the end of each time step.
1480000e7ae3SMatthew Knepley 
14813f9fe445SBarry Smith   Logically Collective on TS
1482000e7ae3SMatthew Knepley 
1483000e7ae3SMatthew Knepley   Input Parameters:
1484000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1485000e7ae3SMatthew Knepley - func - The function
1486000e7ae3SMatthew Knepley 
1487000e7ae3SMatthew Knepley   Calling sequence of func:
1488000e7ae3SMatthew Knepley . func (TS ts);
1489000e7ae3SMatthew Knepley 
1490000e7ae3SMatthew Knepley   Level: intermediate
1491000e7ae3SMatthew Knepley 
1492000e7ae3SMatthew Knepley .keywords: TS, timestep
1493000e7ae3SMatthew Knepley @*/
14947087cfbeSBarry Smith PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1495000e7ae3SMatthew Knepley {
1496000e7ae3SMatthew Knepley   PetscFunctionBegin;
14970700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1498000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1499000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1500000e7ae3SMatthew Knepley }
1501000e7ae3SMatthew Knepley 
1502e74ef692SMatthew Knepley #undef __FUNCT__
15033f2090d5SJed Brown #define __FUNCT__ "TSPostStep"
15043f2090d5SJed Brown /*@C
15053f2090d5SJed Brown   TSPostStep - Runs the user-defined post-step function.
15063f2090d5SJed Brown 
15073f2090d5SJed Brown   Collective on TS
15083f2090d5SJed Brown 
15093f2090d5SJed Brown   Input Parameters:
15103f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
15113f2090d5SJed Brown 
15123f2090d5SJed Brown   Notes:
15133f2090d5SJed Brown   TSPostStep() is typically used within time stepping implementations,
15143f2090d5SJed Brown   so most users would not generally call this routine themselves.
15153f2090d5SJed Brown 
15163f2090d5SJed Brown   Level: developer
15173f2090d5SJed Brown 
15183f2090d5SJed Brown .keywords: TS, timestep
15193f2090d5SJed Brown @*/
15207087cfbeSBarry Smith PetscErrorCode  TSPostStep(TS ts)
15213f2090d5SJed Brown {
15223f2090d5SJed Brown   PetscErrorCode ierr;
15233f2090d5SJed Brown 
15243f2090d5SJed Brown   PetscFunctionBegin;
15250700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
152672ac3e02SJed Brown   if (ts->ops->poststep) {
15273f2090d5SJed Brown     PetscStackPush("TS PostStep function");
15283f2090d5SJed Brown     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
15293f2090d5SJed Brown     PetscStackPop;
153072ac3e02SJed Brown   }
15313f2090d5SJed Brown   PetscFunctionReturn(0);
15323f2090d5SJed Brown }
15333f2090d5SJed Brown 
1534d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1535d763cef2SBarry Smith 
15364a2ae208SSatish Balay #undef __FUNCT__
1537a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet"
1538d763cef2SBarry Smith /*@C
1539a6570f20SBarry Smith    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1540d763cef2SBarry Smith    timestep to display the iteration's  progress.
1541d763cef2SBarry Smith 
15423f9fe445SBarry Smith    Logically Collective on TS
1543d763cef2SBarry Smith 
1544d763cef2SBarry Smith    Input Parameters:
1545d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1546d763cef2SBarry Smith .  func - monitoring routine
1547329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1548b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1549b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1550b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1551d763cef2SBarry Smith 
1552d763cef2SBarry Smith    Calling sequence of func:
1553a7cc72afSBarry Smith $    int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1554d763cef2SBarry Smith 
1555d763cef2SBarry Smith +    ts - the TS context
1556d763cef2SBarry Smith .    steps - iteration number
15571f06c33eSBarry Smith .    time - current time
1558d763cef2SBarry Smith .    x - current iterate
1559d763cef2SBarry Smith -    mctx - [optional] monitoring context
1560d763cef2SBarry Smith 
1561d763cef2SBarry Smith    Notes:
1562d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1563d763cef2SBarry Smith    already has been loaded.
1564d763cef2SBarry Smith 
1565025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each TS object
1566025f1a04SBarry Smith 
1567d763cef2SBarry Smith    Level: intermediate
1568d763cef2SBarry Smith 
1569d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1570d763cef2SBarry Smith 
1571a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel()
1572d763cef2SBarry Smith @*/
1573c2efdce3SBarry Smith PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1574d763cef2SBarry Smith {
1575d763cef2SBarry Smith   PetscFunctionBegin;
15760700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
157717186662SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1578d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1579329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1580d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1581d763cef2SBarry Smith   PetscFunctionReturn(0);
1582d763cef2SBarry Smith }
1583d763cef2SBarry Smith 
15844a2ae208SSatish Balay #undef __FUNCT__
1585a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel"
1586d763cef2SBarry Smith /*@C
1587a6570f20SBarry Smith    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1588d763cef2SBarry Smith 
15893f9fe445SBarry Smith    Logically Collective on TS
1590d763cef2SBarry Smith 
1591d763cef2SBarry Smith    Input Parameters:
1592d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1593d763cef2SBarry Smith 
1594d763cef2SBarry Smith    Notes:
1595d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1596d763cef2SBarry Smith 
1597d763cef2SBarry Smith    Level: intermediate
1598d763cef2SBarry Smith 
1599d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1600d763cef2SBarry Smith 
1601a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet()
1602d763cef2SBarry Smith @*/
16037087cfbeSBarry Smith PetscErrorCode  TSMonitorCancel(TS ts)
1604d763cef2SBarry Smith {
1605d952e501SBarry Smith   PetscErrorCode ierr;
1606d952e501SBarry Smith   PetscInt       i;
1607d952e501SBarry Smith 
1608d763cef2SBarry Smith   PetscFunctionBegin;
16090700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1610d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1611d952e501SBarry Smith     if (ts->mdestroy[i]) {
16123c4aec1bSBarry Smith       ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
1613d952e501SBarry Smith     }
1614d952e501SBarry Smith   }
1615d763cef2SBarry Smith   ts->numbermonitors = 0;
1616d763cef2SBarry Smith   PetscFunctionReturn(0);
1617d763cef2SBarry Smith }
1618d763cef2SBarry Smith 
16194a2ae208SSatish Balay #undef __FUNCT__
1620a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault"
1621d8e5e3e6SSatish Balay /*@
1622a6570f20SBarry Smith    TSMonitorDefault - Sets the Default monitor
16235516499fSSatish Balay 
16245516499fSSatish Balay    Level: intermediate
162541251cbbSSatish Balay 
16265516499fSSatish Balay .keywords: TS, set, monitor
16275516499fSSatish Balay 
162841251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet()
162941251cbbSSatish Balay @*/
1630649052a6SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1631d763cef2SBarry Smith {
1632dfbe8321SBarry Smith   PetscErrorCode ierr;
1633649052a6SBarry Smith   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1634d132466eSBarry Smith 
1635d763cef2SBarry Smith   PetscFunctionBegin;
1636649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1637649052a6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr);
1638649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1639d763cef2SBarry Smith   PetscFunctionReturn(0);
1640d763cef2SBarry Smith }
1641d763cef2SBarry Smith 
16424a2ae208SSatish Balay #undef __FUNCT__
16434a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1644d763cef2SBarry Smith /*@
1645d763cef2SBarry Smith    TSStep - Steps the requested number of timesteps.
1646d763cef2SBarry Smith 
1647d763cef2SBarry Smith    Collective on TS
1648d763cef2SBarry Smith 
1649d763cef2SBarry Smith    Input Parameter:
1650d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1651d763cef2SBarry Smith 
1652d763cef2SBarry Smith    Output Parameters:
1653d763cef2SBarry Smith +  steps - number of iterations until termination
1654142b95e3SSatish Balay -  ptime - time until termination
1655d763cef2SBarry Smith 
1656d763cef2SBarry Smith    Level: beginner
1657d763cef2SBarry Smith 
1658d763cef2SBarry Smith .keywords: TS, timestep, solve
1659d763cef2SBarry Smith 
1660d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSDestroy()
1661d763cef2SBarry Smith @*/
1662193ac0bcSJed Brown PetscErrorCode  TSStep(TS ts)
1663d763cef2SBarry Smith {
1664dfbe8321SBarry Smith   PetscErrorCode ierr;
1665d763cef2SBarry Smith 
1666d763cef2SBarry Smith   PetscFunctionBegin;
16670700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1668277b19d0SLisandro Dalcin 
1669d405a339SMatthew Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
1670d405a339SMatthew Knepley 
1671d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1672193ac0bcSJed Brown   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
1673d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr);
1674d763cef2SBarry Smith   PetscFunctionReturn(0);
1675d763cef2SBarry Smith }
1676d763cef2SBarry Smith 
16774a2ae208SSatish Balay #undef __FUNCT__
16786a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve"
16796a4d4014SLisandro Dalcin /*@
16806a4d4014SLisandro Dalcin    TSSolve - Steps the requested number of timesteps.
16816a4d4014SLisandro Dalcin 
16826a4d4014SLisandro Dalcin    Collective on TS
16836a4d4014SLisandro Dalcin 
16846a4d4014SLisandro Dalcin    Input Parameter:
16856a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
16866a4d4014SLisandro Dalcin -  x - the solution vector, or PETSC_NULL if it was set with TSSetSolution()
16876a4d4014SLisandro Dalcin 
16886a4d4014SLisandro Dalcin    Level: beginner
16896a4d4014SLisandro Dalcin 
16906a4d4014SLisandro Dalcin .keywords: TS, timestep, solve
16916a4d4014SLisandro Dalcin 
16926a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep()
16936a4d4014SLisandro Dalcin @*/
16947087cfbeSBarry Smith PetscErrorCode  TSSolve(TS ts, Vec x)
16956a4d4014SLisandro Dalcin {
1696193ac0bcSJed Brown   PetscInt       i;
16976a4d4014SLisandro Dalcin   PetscErrorCode ierr;
1698f22f69f0SBarry Smith 
16996a4d4014SLisandro Dalcin   PetscFunctionBegin;
17000700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1701193ac0bcSJed Brown   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1702193ac0bcSJed Brown   ierr = TSSetSolution(ts,x); CHKERRQ(ierr);
17036a4d4014SLisandro Dalcin   /* reset time step and iteration counters */
1704193ac0bcSJed Brown   ts->steps = 0;
1705193ac0bcSJed Brown   ts->linear_its = 0;
1706193ac0bcSJed Brown   ts->nonlinear_its = 0;
1707193ac0bcSJed Brown   ts->reason = TS_CONVERGED_ITERATING;
1708193ac0bcSJed Brown   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1709193ac0bcSJed Brown 
1710193ac0bcSJed Brown   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
1711193ac0bcSJed Brown     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
1712193ac0bcSJed Brown   } else {
17136a4d4014SLisandro Dalcin     /* steps the requested number of timesteps. */
1714*3daf2de8SJed Brown     for (i=0; !ts->reason; ) {
1715193ac0bcSJed Brown       ierr = TSPreStep(ts);CHKERRQ(ierr);
1716193ac0bcSJed Brown       ierr = TSStep(ts);CHKERRQ(ierr);
1717193ac0bcSJed Brown       if (ts->reason < 0) {
1718193ac0bcSJed Brown         if (ts->errorifstepfailed) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed");
1719*3daf2de8SJed Brown       } else if (++i >= ts->max_steps) {
1720193ac0bcSJed Brown         ts->reason = TS_CONVERGED_ITS;
1721193ac0bcSJed Brown       } else if (ts->ptime >= ts->max_time) {
1722193ac0bcSJed Brown         ts->reason = TS_CONVERGED_TIME;
1723193ac0bcSJed Brown       }
1724193ac0bcSJed Brown       ierr = TSPostStep(ts);CHKERRQ(ierr);
1725193ac0bcSJed Brown       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1726193ac0bcSJed Brown     }
1727193ac0bcSJed Brown   }
1728193ac0bcSJed Brown   if (!PetscPreLoadingOn) {
1729193ac0bcSJed Brown     ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr);
1730193ac0bcSJed Brown   }
17316a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
17326a4d4014SLisandro Dalcin }
17336a4d4014SLisandro Dalcin 
17346a4d4014SLisandro Dalcin #undef __FUNCT__
17354a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1736d763cef2SBarry Smith /*
1737d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1738d763cef2SBarry Smith */
1739a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1740d763cef2SBarry Smith {
17416849ba73SBarry Smith   PetscErrorCode ierr;
1742a7cc72afSBarry Smith   PetscInt       i,n = ts->numbermonitors;
1743d763cef2SBarry Smith 
1744d763cef2SBarry Smith   PetscFunctionBegin;
1745d763cef2SBarry Smith   for (i=0; i<n; i++) {
1746142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1747d763cef2SBarry Smith   }
1748d763cef2SBarry Smith   PetscFunctionReturn(0);
1749d763cef2SBarry Smith }
1750d763cef2SBarry Smith 
1751d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1752d763cef2SBarry Smith 
17534a2ae208SSatish Balay #undef __FUNCT__
1754a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate"
1755d763cef2SBarry Smith /*@C
1756a6570f20SBarry Smith    TSMonitorLGCreate - Creates a line graph context for use with
1757d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1758d763cef2SBarry Smith 
1759d763cef2SBarry Smith    Collective on TS
1760d763cef2SBarry Smith 
1761d763cef2SBarry Smith    Input Parameters:
1762d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1763d763cef2SBarry Smith .  label - the title to put in the title bar
17647c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1765d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1766d763cef2SBarry Smith 
1767d763cef2SBarry Smith    Output Parameter:
1768d763cef2SBarry Smith .  draw - the drawing context
1769d763cef2SBarry Smith 
1770d763cef2SBarry Smith    Options Database Key:
1771a6570f20SBarry Smith .  -ts_monitor_draw - automatically sets line graph monitor
1772d763cef2SBarry Smith 
1773d763cef2SBarry Smith    Notes:
1774a6570f20SBarry Smith    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1775d763cef2SBarry Smith 
1776d763cef2SBarry Smith    Level: intermediate
1777d763cef2SBarry Smith 
17787c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1779d763cef2SBarry Smith 
1780a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet()
17817c922b88SBarry Smith 
1782d763cef2SBarry Smith @*/
17837087cfbeSBarry Smith PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1784d763cef2SBarry Smith {
1785b0a32e0cSBarry Smith   PetscDraw      win;
1786dfbe8321SBarry Smith   PetscErrorCode ierr;
1787d763cef2SBarry Smith 
1788d763cef2SBarry Smith   PetscFunctionBegin;
1789b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1790b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1791b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1792b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1793d763cef2SBarry Smith 
179452e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1795d763cef2SBarry Smith   PetscFunctionReturn(0);
1796d763cef2SBarry Smith }
1797d763cef2SBarry Smith 
17984a2ae208SSatish Balay #undef __FUNCT__
1799a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG"
1800a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1801d763cef2SBarry Smith {
1802b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
180387828ca2SBarry Smith   PetscReal      x,y = ptime;
1804dfbe8321SBarry Smith   PetscErrorCode ierr;
1805d763cef2SBarry Smith 
1806d763cef2SBarry Smith   PetscFunctionBegin;
18077c922b88SBarry Smith   if (!monctx) {
18087c922b88SBarry Smith     MPI_Comm    comm;
1809b0a32e0cSBarry Smith     PetscViewer viewer;
18107c922b88SBarry Smith 
18117c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1812b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1813b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
18147c922b88SBarry Smith   }
18157c922b88SBarry Smith 
1816b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
181787828ca2SBarry Smith   x = (PetscReal)n;
1818b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1819d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1820b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1821d763cef2SBarry Smith   }
1822d763cef2SBarry Smith   PetscFunctionReturn(0);
1823d763cef2SBarry Smith }
1824d763cef2SBarry Smith 
18254a2ae208SSatish Balay #undef __FUNCT__
1826a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy"
1827d763cef2SBarry Smith /*@C
1828a6570f20SBarry Smith    TSMonitorLGDestroy - Destroys a line graph context that was created
1829a6570f20SBarry Smith    with TSMonitorLGCreate().
1830d763cef2SBarry Smith 
1831b0a32e0cSBarry Smith    Collective on PetscDrawLG
1832d763cef2SBarry Smith 
1833d763cef2SBarry Smith    Input Parameter:
1834d763cef2SBarry Smith .  draw - the drawing context
1835d763cef2SBarry Smith 
1836d763cef2SBarry Smith    Level: intermediate
1837d763cef2SBarry Smith 
1838d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1839d763cef2SBarry Smith 
1840a6570f20SBarry Smith .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1841d763cef2SBarry Smith @*/
18426bf464f9SBarry Smith PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
1843d763cef2SBarry Smith {
1844b0a32e0cSBarry Smith   PetscDraw      draw;
1845dfbe8321SBarry Smith   PetscErrorCode ierr;
1846d763cef2SBarry Smith 
1847d763cef2SBarry Smith   PetscFunctionBegin;
18486bf464f9SBarry Smith   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
18496bf464f9SBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
1850b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1851d763cef2SBarry Smith   PetscFunctionReturn(0);
1852d763cef2SBarry Smith }
1853d763cef2SBarry Smith 
18544a2ae208SSatish Balay #undef __FUNCT__
18554a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
1856d763cef2SBarry Smith /*@
1857d763cef2SBarry Smith    TSGetTime - Gets the current time.
1858d763cef2SBarry Smith 
1859d763cef2SBarry Smith    Not Collective
1860d763cef2SBarry Smith 
1861d763cef2SBarry Smith    Input Parameter:
1862d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1863d763cef2SBarry Smith 
1864d763cef2SBarry Smith    Output Parameter:
1865d763cef2SBarry Smith .  t  - the current time
1866d763cef2SBarry Smith 
1867d763cef2SBarry Smith    Level: beginner
1868d763cef2SBarry Smith 
1869d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1870d763cef2SBarry Smith 
1871d763cef2SBarry Smith .keywords: TS, get, time
1872d763cef2SBarry Smith @*/
18737087cfbeSBarry Smith PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
1874d763cef2SBarry Smith {
1875d763cef2SBarry Smith   PetscFunctionBegin;
18760700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
18774482741eSBarry Smith   PetscValidDoublePointer(t,2);
1878d763cef2SBarry Smith   *t = ts->ptime;
1879d763cef2SBarry Smith   PetscFunctionReturn(0);
1880d763cef2SBarry Smith }
1881d763cef2SBarry Smith 
18824a2ae208SSatish Balay #undef __FUNCT__
18836a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime"
18846a4d4014SLisandro Dalcin /*@
18856a4d4014SLisandro Dalcin    TSSetTime - Allows one to reset the time.
18866a4d4014SLisandro Dalcin 
18873f9fe445SBarry Smith    Logically Collective on TS
18886a4d4014SLisandro Dalcin 
18896a4d4014SLisandro Dalcin    Input Parameters:
18906a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
18916a4d4014SLisandro Dalcin -  time - the time
18926a4d4014SLisandro Dalcin 
18936a4d4014SLisandro Dalcin    Level: intermediate
18946a4d4014SLisandro Dalcin 
18956a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration()
18966a4d4014SLisandro Dalcin 
18976a4d4014SLisandro Dalcin .keywords: TS, set, time
18986a4d4014SLisandro Dalcin @*/
18997087cfbeSBarry Smith PetscErrorCode  TSSetTime(TS ts, PetscReal t)
19006a4d4014SLisandro Dalcin {
19016a4d4014SLisandro Dalcin   PetscFunctionBegin;
19020700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1903c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,t,2);
19046a4d4014SLisandro Dalcin   ts->ptime = t;
19056a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
19066a4d4014SLisandro Dalcin }
19076a4d4014SLisandro Dalcin 
19086a4d4014SLisandro Dalcin #undef __FUNCT__
19094a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
1910d763cef2SBarry Smith /*@C
1911d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
1912d763cef2SBarry Smith    TS options in the database.
1913d763cef2SBarry Smith 
19143f9fe445SBarry Smith    Logically Collective on TS
1915d763cef2SBarry Smith 
1916d763cef2SBarry Smith    Input Parameter:
1917d763cef2SBarry Smith +  ts     - The TS context
1918d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1919d763cef2SBarry Smith 
1920d763cef2SBarry Smith    Notes:
1921d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1922d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1923d763cef2SBarry Smith    hyphen.
1924d763cef2SBarry Smith 
1925d763cef2SBarry Smith    Level: advanced
1926d763cef2SBarry Smith 
1927d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
1928d763cef2SBarry Smith 
1929d763cef2SBarry Smith .seealso: TSSetFromOptions()
1930d763cef2SBarry Smith 
1931d763cef2SBarry Smith @*/
19327087cfbeSBarry Smith PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
1933d763cef2SBarry Smith {
1934dfbe8321SBarry Smith   PetscErrorCode ierr;
1935089b2837SJed Brown   SNES           snes;
1936d763cef2SBarry Smith 
1937d763cef2SBarry Smith   PetscFunctionBegin;
19380700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1939d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1940089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1941089b2837SJed Brown   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
1942d763cef2SBarry Smith   PetscFunctionReturn(0);
1943d763cef2SBarry Smith }
1944d763cef2SBarry Smith 
1945d763cef2SBarry Smith 
19464a2ae208SSatish Balay #undef __FUNCT__
19474a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
1948d763cef2SBarry Smith /*@C
1949d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
1950d763cef2SBarry Smith    TS options in the database.
1951d763cef2SBarry Smith 
19523f9fe445SBarry Smith    Logically Collective on TS
1953d763cef2SBarry Smith 
1954d763cef2SBarry Smith    Input Parameter:
1955d763cef2SBarry Smith +  ts     - The TS context
1956d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
1957d763cef2SBarry Smith 
1958d763cef2SBarry Smith    Notes:
1959d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
1960d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
1961d763cef2SBarry Smith    hyphen.
1962d763cef2SBarry Smith 
1963d763cef2SBarry Smith    Level: advanced
1964d763cef2SBarry Smith 
1965d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
1966d763cef2SBarry Smith 
1967d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
1968d763cef2SBarry Smith 
1969d763cef2SBarry Smith @*/
19707087cfbeSBarry Smith PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
1971d763cef2SBarry Smith {
1972dfbe8321SBarry Smith   PetscErrorCode ierr;
1973089b2837SJed Brown   SNES           snes;
1974d763cef2SBarry Smith 
1975d763cef2SBarry Smith   PetscFunctionBegin;
19760700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1977d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
1978089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1979089b2837SJed Brown   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
1980d763cef2SBarry Smith   PetscFunctionReturn(0);
1981d763cef2SBarry Smith }
1982d763cef2SBarry Smith 
19834a2ae208SSatish Balay #undef __FUNCT__
19844a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
1985d763cef2SBarry Smith /*@C
1986d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
1987d763cef2SBarry Smith    TS options in the database.
1988d763cef2SBarry Smith 
1989d763cef2SBarry Smith    Not Collective
1990d763cef2SBarry Smith 
1991d763cef2SBarry Smith    Input Parameter:
1992d763cef2SBarry Smith .  ts - The TS context
1993d763cef2SBarry Smith 
1994d763cef2SBarry Smith    Output Parameter:
1995d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
1996d763cef2SBarry Smith 
1997d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
1998d763cef2SBarry Smith    sufficient length to hold the prefix.
1999d763cef2SBarry Smith 
2000d763cef2SBarry Smith    Level: intermediate
2001d763cef2SBarry Smith 
2002d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
2003d763cef2SBarry Smith 
2004d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
2005d763cef2SBarry Smith @*/
20067087cfbeSBarry Smith PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2007d763cef2SBarry Smith {
2008dfbe8321SBarry Smith   PetscErrorCode ierr;
2009d763cef2SBarry Smith 
2010d763cef2SBarry Smith   PetscFunctionBegin;
20110700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
20124482741eSBarry Smith   PetscValidPointer(prefix,2);
2013d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2014d763cef2SBarry Smith   PetscFunctionReturn(0);
2015d763cef2SBarry Smith }
2016d763cef2SBarry Smith 
20174a2ae208SSatish Balay #undef __FUNCT__
20184a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
2019d763cef2SBarry Smith /*@C
2020d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2021d763cef2SBarry Smith 
2022d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
2023d763cef2SBarry Smith 
2024d763cef2SBarry Smith    Input Parameter:
2025d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
2026d763cef2SBarry Smith 
2027d763cef2SBarry Smith    Output Parameters:
2028d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
2029d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
2030089b2837SJed Brown .  func - Function to compute the Jacobian of the RHS
2031d763cef2SBarry Smith -  ctx - User-defined context for Jacobian evaluation routine
2032d763cef2SBarry Smith 
2033d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2034d763cef2SBarry Smith 
2035d763cef2SBarry Smith    Level: intermediate
2036d763cef2SBarry Smith 
203726d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2038d763cef2SBarry Smith 
2039d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
2040d763cef2SBarry Smith @*/
2041089b2837SJed Brown PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2042d763cef2SBarry Smith {
2043089b2837SJed Brown   PetscErrorCode ierr;
2044089b2837SJed Brown   SNES           snes;
2045089b2837SJed Brown 
2046d763cef2SBarry Smith   PetscFunctionBegin;
2047089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2048089b2837SJed Brown   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2049089b2837SJed Brown   if (func) *func = ts->ops->rhsjacobian;
205026d46c62SHong Zhang   if (ctx) *ctx = ts->jacP;
2051d763cef2SBarry Smith   PetscFunctionReturn(0);
2052d763cef2SBarry Smith }
2053d763cef2SBarry Smith 
20541713a123SBarry Smith #undef __FUNCT__
20552eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian"
20562eca1d9cSJed Brown /*@C
20572eca1d9cSJed Brown    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
20582eca1d9cSJed Brown 
20592eca1d9cSJed Brown    Not Collective, but parallel objects are returned if TS is parallel
20602eca1d9cSJed Brown 
20612eca1d9cSJed Brown    Input Parameter:
20622eca1d9cSJed Brown .  ts  - The TS context obtained from TSCreate()
20632eca1d9cSJed Brown 
20642eca1d9cSJed Brown    Output Parameters:
20652eca1d9cSJed Brown +  A   - The Jacobian of F(t,U,U_t)
20662eca1d9cSJed Brown .  B   - The preconditioner matrix, often the same as A
20672eca1d9cSJed Brown .  f   - The function to compute the matrices
20682eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine
20692eca1d9cSJed Brown 
20702eca1d9cSJed Brown    Notes: You can pass in PETSC_NULL for any return argument you do not need.
20712eca1d9cSJed Brown 
20722eca1d9cSJed Brown    Level: advanced
20732eca1d9cSJed Brown 
20742eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
20752eca1d9cSJed Brown 
20762eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian
20772eca1d9cSJed Brown @*/
20787087cfbeSBarry Smith PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
20792eca1d9cSJed Brown {
2080089b2837SJed Brown   PetscErrorCode ierr;
2081089b2837SJed Brown   SNES           snes;
2082089b2837SJed Brown 
20832eca1d9cSJed Brown   PetscFunctionBegin;
2084089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2085089b2837SJed Brown   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
20862eca1d9cSJed Brown   if (f) *f = ts->ops->ijacobian;
20872eca1d9cSJed Brown   if (ctx) *ctx = ts->jacP;
20882eca1d9cSJed Brown   PetscFunctionReturn(0);
20892eca1d9cSJed Brown }
20902eca1d9cSJed Brown 
20912eca1d9cSJed Brown #undef __FUNCT__
2092a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution"
20931713a123SBarry Smith /*@C
2094a6570f20SBarry Smith    TSMonitorSolution - Monitors progress of the TS solvers by calling
20951713a123SBarry Smith    VecView() for the solution at each timestep
20961713a123SBarry Smith 
20971713a123SBarry Smith    Collective on TS
20981713a123SBarry Smith 
20991713a123SBarry Smith    Input Parameters:
21001713a123SBarry Smith +  ts - the TS context
21011713a123SBarry Smith .  step - current time-step
2102142b95e3SSatish Balay .  ptime - current time
21031713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
21041713a123SBarry Smith 
21051713a123SBarry Smith    Level: intermediate
21061713a123SBarry Smith 
21071713a123SBarry Smith .keywords: TS,  vector, monitor, view
21081713a123SBarry Smith 
2109a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
21101713a123SBarry Smith @*/
21117087cfbeSBarry Smith PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
21121713a123SBarry Smith {
2113dfbe8321SBarry Smith   PetscErrorCode ierr;
21141713a123SBarry Smith   PetscViewer    viewer = (PetscViewer) dummy;
21151713a123SBarry Smith 
21161713a123SBarry Smith   PetscFunctionBegin;
2117a34d58ebSBarry Smith   if (!dummy) {
21187adad957SLisandro Dalcin     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
21191713a123SBarry Smith   }
21201713a123SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
21211713a123SBarry Smith   PetscFunctionReturn(0);
21221713a123SBarry Smith }
21231713a123SBarry Smith 
21241713a123SBarry Smith 
21256c699258SBarry Smith #undef __FUNCT__
21266c699258SBarry Smith #define __FUNCT__ "TSSetDM"
21276c699258SBarry Smith /*@
21286c699258SBarry Smith    TSSetDM - Sets the DM that may be used by some preconditioners
21296c699258SBarry Smith 
21303f9fe445SBarry Smith    Logically Collective on TS and DM
21316c699258SBarry Smith 
21326c699258SBarry Smith    Input Parameters:
21336c699258SBarry Smith +  ts - the preconditioner context
21346c699258SBarry Smith -  dm - the dm
21356c699258SBarry Smith 
21366c699258SBarry Smith    Level: intermediate
21376c699258SBarry Smith 
21386c699258SBarry Smith 
21396c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
21406c699258SBarry Smith @*/
21417087cfbeSBarry Smith PetscErrorCode  TSSetDM(TS ts,DM dm)
21426c699258SBarry Smith {
21436c699258SBarry Smith   PetscErrorCode ierr;
2144089b2837SJed Brown   SNES           snes;
21456c699258SBarry Smith 
21466c699258SBarry Smith   PetscFunctionBegin;
21470700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
214870663e4aSLisandro Dalcin   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
21496bf464f9SBarry Smith   ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
21506c699258SBarry Smith   ts->dm = dm;
2151089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2152089b2837SJed Brown   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
21536c699258SBarry Smith   PetscFunctionReturn(0);
21546c699258SBarry Smith }
21556c699258SBarry Smith 
21566c699258SBarry Smith #undef __FUNCT__
21576c699258SBarry Smith #define __FUNCT__ "TSGetDM"
21586c699258SBarry Smith /*@
21596c699258SBarry Smith    TSGetDM - Gets the DM that may be used by some preconditioners
21606c699258SBarry Smith 
21613f9fe445SBarry Smith    Not Collective
21626c699258SBarry Smith 
21636c699258SBarry Smith    Input Parameter:
21646c699258SBarry Smith . ts - the preconditioner context
21656c699258SBarry Smith 
21666c699258SBarry Smith    Output Parameter:
21676c699258SBarry Smith .  dm - the dm
21686c699258SBarry Smith 
21696c699258SBarry Smith    Level: intermediate
21706c699258SBarry Smith 
21716c699258SBarry Smith 
21726c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
21736c699258SBarry Smith @*/
21747087cfbeSBarry Smith PetscErrorCode  TSGetDM(TS ts,DM *dm)
21756c699258SBarry Smith {
21766c699258SBarry Smith   PetscFunctionBegin;
21770700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
21786c699258SBarry Smith   *dm = ts->dm;
21796c699258SBarry Smith   PetscFunctionReturn(0);
21806c699258SBarry Smith }
21811713a123SBarry Smith 
21820f5c6efeSJed Brown #undef __FUNCT__
21830f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction"
21840f5c6efeSJed Brown /*@
21850f5c6efeSJed Brown    SNESTSFormFunction - Function to evaluate nonlinear residual
21860f5c6efeSJed Brown 
21873f9fe445SBarry Smith    Logically Collective on SNES
21880f5c6efeSJed Brown 
21890f5c6efeSJed Brown    Input Parameter:
2190d42a1c89SJed Brown + snes - nonlinear solver
21910f5c6efeSJed Brown . X - the current state at which to evaluate the residual
2192d42a1c89SJed Brown - ctx - user context, must be a TS
21930f5c6efeSJed Brown 
21940f5c6efeSJed Brown    Output Parameter:
21950f5c6efeSJed Brown . F - the nonlinear residual
21960f5c6efeSJed Brown 
21970f5c6efeSJed Brown    Notes:
21980f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
21990f5c6efeSJed Brown    It is most frequently passed to MatFDColoringSetFunction().
22000f5c6efeSJed Brown 
22010f5c6efeSJed Brown    Level: advanced
22020f5c6efeSJed Brown 
22030f5c6efeSJed Brown .seealso: SNESSetFunction(), MatFDColoringSetFunction()
22040f5c6efeSJed Brown @*/
22057087cfbeSBarry Smith PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
22060f5c6efeSJed Brown {
22070f5c6efeSJed Brown   TS ts = (TS)ctx;
22080f5c6efeSJed Brown   PetscErrorCode ierr;
22090f5c6efeSJed Brown 
22100f5c6efeSJed Brown   PetscFunctionBegin;
22110f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
22120f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
22130f5c6efeSJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
22140f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
22150f5c6efeSJed Brown   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
22160f5c6efeSJed Brown   PetscFunctionReturn(0);
22170f5c6efeSJed Brown }
22180f5c6efeSJed Brown 
22190f5c6efeSJed Brown #undef __FUNCT__
22200f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian"
22210f5c6efeSJed Brown /*@
22220f5c6efeSJed Brown    SNESTSFormJacobian - Function to evaluate the Jacobian
22230f5c6efeSJed Brown 
22240f5c6efeSJed Brown    Collective on SNES
22250f5c6efeSJed Brown 
22260f5c6efeSJed Brown    Input Parameter:
22270f5c6efeSJed Brown + snes - nonlinear solver
22280f5c6efeSJed Brown . X - the current state at which to evaluate the residual
22290f5c6efeSJed Brown - ctx - user context, must be a TS
22300f5c6efeSJed Brown 
22310f5c6efeSJed Brown    Output Parameter:
22320f5c6efeSJed Brown + A - the Jacobian
22330f5c6efeSJed Brown . B - the preconditioning matrix (may be the same as A)
22340f5c6efeSJed Brown - flag - indicates any structure change in the matrix
22350f5c6efeSJed Brown 
22360f5c6efeSJed Brown    Notes:
22370f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
22380f5c6efeSJed Brown 
22390f5c6efeSJed Brown    Level: developer
22400f5c6efeSJed Brown 
22410f5c6efeSJed Brown .seealso: SNESSetJacobian()
22420f5c6efeSJed Brown @*/
22437087cfbeSBarry Smith PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
22440f5c6efeSJed Brown {
22450f5c6efeSJed Brown   TS ts = (TS)ctx;
22460f5c6efeSJed Brown   PetscErrorCode ierr;
22470f5c6efeSJed Brown 
22480f5c6efeSJed Brown   PetscFunctionBegin;
22490f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
22500f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
22510f5c6efeSJed Brown   PetscValidPointer(A,3);
22520f5c6efeSJed Brown   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
22530f5c6efeSJed Brown   PetscValidPointer(B,4);
22540f5c6efeSJed Brown   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
22550f5c6efeSJed Brown   PetscValidPointer(flag,5);
22560f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
22570f5c6efeSJed Brown   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
22580f5c6efeSJed Brown   PetscFunctionReturn(0);
22590f5c6efeSJed Brown }
2260325fc9f4SBarry Smith 
2261325fc9f4SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
2262c6db04a5SJed Brown #include <mex.h>
2263325fc9f4SBarry Smith 
2264325fc9f4SBarry Smith typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
2265325fc9f4SBarry Smith 
2266325fc9f4SBarry Smith #undef __FUNCT__
2267325fc9f4SBarry Smith #define __FUNCT__ "TSComputeFunction_Matlab"
2268325fc9f4SBarry Smith /*
2269325fc9f4SBarry Smith    TSComputeFunction_Matlab - Calls the function that has been set with
2270325fc9f4SBarry Smith                          TSSetFunctionMatlab().
2271325fc9f4SBarry Smith 
2272325fc9f4SBarry Smith    Collective on TS
2273325fc9f4SBarry Smith 
2274325fc9f4SBarry Smith    Input Parameters:
2275325fc9f4SBarry Smith +  snes - the TS context
2276325fc9f4SBarry Smith -  x - input vector
2277325fc9f4SBarry Smith 
2278325fc9f4SBarry Smith    Output Parameter:
2279325fc9f4SBarry Smith .  y - function vector, as set by TSSetFunction()
2280325fc9f4SBarry Smith 
2281325fc9f4SBarry Smith    Notes:
2282325fc9f4SBarry Smith    TSComputeFunction() is typically used within nonlinear solvers
2283325fc9f4SBarry Smith    implementations, so most users would not generally call this routine
2284325fc9f4SBarry Smith    themselves.
2285325fc9f4SBarry Smith 
2286325fc9f4SBarry Smith    Level: developer
2287325fc9f4SBarry Smith 
2288325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
2289325fc9f4SBarry Smith 
2290325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
2291325fc9f4SBarry Smith */
22927087cfbeSBarry Smith PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
2293325fc9f4SBarry Smith {
2294325fc9f4SBarry Smith   PetscErrorCode   ierr;
2295325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2296325fc9f4SBarry Smith   int              nlhs = 1,nrhs = 7;
2297325fc9f4SBarry Smith   mxArray	   *plhs[1],*prhs[7];
2298325fc9f4SBarry Smith   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
2299325fc9f4SBarry Smith 
2300325fc9f4SBarry Smith   PetscFunctionBegin;
2301325fc9f4SBarry Smith   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2302325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2303325fc9f4SBarry Smith   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
2304325fc9f4SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
2305325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,x,3);
2306325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,y,5);
2307325fc9f4SBarry Smith 
2308325fc9f4SBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2309325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2310880f3077SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
2311325fc9f4SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
2312325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
2313325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar(time);
2314325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
2315325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2316325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)ly);
2317325fc9f4SBarry Smith   prhs[5] =  mxCreateString(sctx->funcname);
2318325fc9f4SBarry Smith   prhs[6] =  sctx->ctx;
2319325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
2320325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2321325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
2322325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
2323325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
2324325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
2325325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
2326325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
2327325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
2328325fc9f4SBarry Smith   PetscFunctionReturn(0);
2329325fc9f4SBarry Smith }
2330325fc9f4SBarry Smith 
2331325fc9f4SBarry Smith 
2332325fc9f4SBarry Smith #undef __FUNCT__
2333325fc9f4SBarry Smith #define __FUNCT__ "TSSetFunctionMatlab"
2334325fc9f4SBarry Smith /*
2335325fc9f4SBarry Smith    TSSetFunctionMatlab - Sets the function evaluation routine and function
2336325fc9f4SBarry Smith    vector for use by the TS routines in solving ODEs
2337e3c5b3baSBarry Smith    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
2338325fc9f4SBarry Smith 
2339325fc9f4SBarry Smith    Logically Collective on TS
2340325fc9f4SBarry Smith 
2341325fc9f4SBarry Smith    Input Parameters:
2342325fc9f4SBarry Smith +  ts - the TS context
2343325fc9f4SBarry Smith -  func - function evaluation routine
2344325fc9f4SBarry Smith 
2345325fc9f4SBarry Smith    Calling sequence of func:
2346325fc9f4SBarry Smith $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
2347325fc9f4SBarry Smith 
2348325fc9f4SBarry Smith    Level: beginner
2349325fc9f4SBarry Smith 
2350325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
2351325fc9f4SBarry Smith 
2352325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2353325fc9f4SBarry Smith */
23547087cfbeSBarry Smith PetscErrorCode  TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx)
2355325fc9f4SBarry Smith {
2356325fc9f4SBarry Smith   PetscErrorCode  ierr;
2357325fc9f4SBarry Smith   TSMatlabContext *sctx;
2358325fc9f4SBarry Smith 
2359325fc9f4SBarry Smith   PetscFunctionBegin;
2360325fc9f4SBarry Smith   /* currently sctx is memory bleed */
2361325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2362325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2363325fc9f4SBarry Smith   /*
2364325fc9f4SBarry Smith      This should work, but it doesn't
2365325fc9f4SBarry Smith   sctx->ctx = ctx;
2366325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
2367325fc9f4SBarry Smith   */
2368325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
2369325fc9f4SBarry Smith   ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
2370325fc9f4SBarry Smith   PetscFunctionReturn(0);
2371325fc9f4SBarry Smith }
2372325fc9f4SBarry Smith 
2373325fc9f4SBarry Smith #undef __FUNCT__
2374325fc9f4SBarry Smith #define __FUNCT__ "TSComputeJacobian_Matlab"
2375325fc9f4SBarry Smith /*
2376325fc9f4SBarry Smith    TSComputeJacobian_Matlab - Calls the function that has been set with
2377325fc9f4SBarry Smith                          TSSetJacobianMatlab().
2378325fc9f4SBarry Smith 
2379325fc9f4SBarry Smith    Collective on TS
2380325fc9f4SBarry Smith 
2381325fc9f4SBarry Smith    Input Parameters:
2382325fc9f4SBarry Smith +  snes - the TS context
2383325fc9f4SBarry Smith .  x - input vector
2384325fc9f4SBarry Smith .  A, B - the matrices
2385325fc9f4SBarry Smith -  ctx - user context
2386325fc9f4SBarry Smith 
2387325fc9f4SBarry Smith    Output Parameter:
2388325fc9f4SBarry Smith .  flag - structure of the matrix
2389325fc9f4SBarry Smith 
2390325fc9f4SBarry Smith    Level: developer
2391325fc9f4SBarry Smith 
2392325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
2393325fc9f4SBarry Smith 
2394325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
2395325fc9f4SBarry Smith @*/
23967087cfbeSBarry Smith PetscErrorCode  TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
2397325fc9f4SBarry Smith {
2398325fc9f4SBarry Smith   PetscErrorCode  ierr;
2399325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2400325fc9f4SBarry Smith   int             nlhs = 2,nrhs = 9;
2401325fc9f4SBarry Smith   mxArray	  *plhs[2],*prhs[9];
2402325fc9f4SBarry Smith   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
2403325fc9f4SBarry Smith 
2404325fc9f4SBarry Smith   PetscFunctionBegin;
2405325fc9f4SBarry Smith   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2406325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2407325fc9f4SBarry Smith 
2408325fc9f4SBarry Smith   /* call Matlab function in ctx with arguments x and y */
2409325fc9f4SBarry Smith 
2410325fc9f4SBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2411325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2412325fc9f4SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
2413325fc9f4SBarry Smith   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
2414325fc9f4SBarry Smith   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
2415325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
2416325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)time);
2417325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
2418325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2419325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)shift);
2420325fc9f4SBarry Smith   prhs[5] =  mxCreateDoubleScalar((double)lA);
2421325fc9f4SBarry Smith   prhs[6] =  mxCreateDoubleScalar((double)lB);
2422325fc9f4SBarry Smith   prhs[7] =  mxCreateString(sctx->funcname);
2423325fc9f4SBarry Smith   prhs[8] =  sctx->ctx;
2424325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
2425325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2426325fc9f4SBarry Smith   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
2427325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
2428325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
2429325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
2430325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
2431325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
2432325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
2433325fc9f4SBarry Smith   mxDestroyArray(prhs[6]);
2434325fc9f4SBarry Smith   mxDestroyArray(prhs[7]);
2435325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
2436325fc9f4SBarry Smith   mxDestroyArray(plhs[1]);
2437325fc9f4SBarry Smith   PetscFunctionReturn(0);
2438325fc9f4SBarry Smith }
2439325fc9f4SBarry Smith 
2440325fc9f4SBarry Smith 
2441325fc9f4SBarry Smith #undef __FUNCT__
2442325fc9f4SBarry Smith #define __FUNCT__ "TSSetJacobianMatlab"
2443325fc9f4SBarry Smith /*
2444325fc9f4SBarry Smith    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
2445e3c5b3baSBarry Smith    vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function
2446325fc9f4SBarry Smith 
2447325fc9f4SBarry Smith    Logically Collective on TS
2448325fc9f4SBarry Smith 
2449325fc9f4SBarry Smith    Input Parameters:
2450325fc9f4SBarry Smith +  snes - the TS context
2451325fc9f4SBarry Smith .  A,B - Jacobian matrices
2452325fc9f4SBarry Smith .  func - function evaluation routine
2453325fc9f4SBarry Smith -  ctx - user context
2454325fc9f4SBarry Smith 
2455325fc9f4SBarry Smith    Calling sequence of func:
2456325fc9f4SBarry Smith $    flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
2457325fc9f4SBarry Smith 
2458325fc9f4SBarry Smith 
2459325fc9f4SBarry Smith    Level: developer
2460325fc9f4SBarry Smith 
2461325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
2462325fc9f4SBarry Smith 
2463325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2464325fc9f4SBarry Smith */
24657087cfbeSBarry Smith PetscErrorCode  TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx)
2466325fc9f4SBarry Smith {
2467325fc9f4SBarry Smith   PetscErrorCode    ierr;
2468325fc9f4SBarry Smith   TSMatlabContext *sctx;
2469325fc9f4SBarry Smith 
2470325fc9f4SBarry Smith   PetscFunctionBegin;
2471325fc9f4SBarry Smith   /* currently sctx is memory bleed */
2472325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2473325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2474325fc9f4SBarry Smith   /*
2475325fc9f4SBarry Smith      This should work, but it doesn't
2476325fc9f4SBarry Smith   sctx->ctx = ctx;
2477325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
2478325fc9f4SBarry Smith   */
2479325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
2480325fc9f4SBarry Smith   ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
2481325fc9f4SBarry Smith   PetscFunctionReturn(0);
2482325fc9f4SBarry Smith }
2483325fc9f4SBarry Smith 
2484b5b1a830SBarry Smith #undef __FUNCT__
2485b5b1a830SBarry Smith #define __FUNCT__ "TSMonitor_Matlab"
2486b5b1a830SBarry Smith /*
2487b5b1a830SBarry Smith    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
2488b5b1a830SBarry Smith 
2489b5b1a830SBarry Smith    Collective on TS
2490b5b1a830SBarry Smith 
2491b5b1a830SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
2492b5b1a830SBarry Smith @*/
24937087cfbeSBarry Smith PetscErrorCode  TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx)
2494b5b1a830SBarry Smith {
2495b5b1a830SBarry Smith   PetscErrorCode  ierr;
2496b5b1a830SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2497a530c242SBarry Smith   int             nlhs = 1,nrhs = 6;
2498b5b1a830SBarry Smith   mxArray	  *plhs[1],*prhs[6];
2499b5b1a830SBarry Smith   long long int   lx = 0,ls = 0;
2500b5b1a830SBarry Smith 
2501b5b1a830SBarry Smith   PetscFunctionBegin;
2502b5b1a830SBarry Smith   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2503b5b1a830SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
2504b5b1a830SBarry Smith 
2505b5b1a830SBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2506b5b1a830SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2507b5b1a830SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
2508b5b1a830SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)it);
2509b5b1a830SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)time);
2510b5b1a830SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lx);
2511b5b1a830SBarry Smith   prhs[4] =  mxCreateString(sctx->funcname);
2512b5b1a830SBarry Smith   prhs[5] =  sctx->ctx;
2513b5b1a830SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
2514b5b1a830SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2515b5b1a830SBarry Smith   mxDestroyArray(prhs[0]);
2516b5b1a830SBarry Smith   mxDestroyArray(prhs[1]);
2517b5b1a830SBarry Smith   mxDestroyArray(prhs[2]);
2518b5b1a830SBarry Smith   mxDestroyArray(prhs[3]);
2519b5b1a830SBarry Smith   mxDestroyArray(prhs[4]);
2520b5b1a830SBarry Smith   mxDestroyArray(plhs[0]);
2521b5b1a830SBarry Smith   PetscFunctionReturn(0);
2522b5b1a830SBarry Smith }
2523b5b1a830SBarry Smith 
2524b5b1a830SBarry Smith 
2525b5b1a830SBarry Smith #undef __FUNCT__
2526b5b1a830SBarry Smith #define __FUNCT__ "TSMonitorSetMatlab"
2527b5b1a830SBarry Smith /*
2528b5b1a830SBarry Smith    TSMonitorSetMatlab - Sets the monitor function from Matlab
2529b5b1a830SBarry Smith 
2530b5b1a830SBarry Smith    Level: developer
2531b5b1a830SBarry Smith 
2532b5b1a830SBarry Smith .keywords: TS, nonlinear, set, function
2533b5b1a830SBarry Smith 
2534b5b1a830SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2535b5b1a830SBarry Smith */
25367087cfbeSBarry Smith PetscErrorCode  TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx)
2537b5b1a830SBarry Smith {
2538b5b1a830SBarry Smith   PetscErrorCode    ierr;
2539b5b1a830SBarry Smith   TSMatlabContext *sctx;
2540b5b1a830SBarry Smith 
2541b5b1a830SBarry Smith   PetscFunctionBegin;
2542b5b1a830SBarry Smith   /* currently sctx is memory bleed */
2543b5b1a830SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2544b5b1a830SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2545b5b1a830SBarry Smith   /*
2546b5b1a830SBarry Smith      This should work, but it doesn't
2547b5b1a830SBarry Smith   sctx->ctx = ctx;
2548b5b1a830SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
2549b5b1a830SBarry Smith   */
2550b5b1a830SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
2551b5b1a830SBarry Smith   ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
2552b5b1a830SBarry Smith   PetscFunctionReturn(0);
2553b5b1a830SBarry Smith }
2554b5b1a830SBarry Smith 
2555325fc9f4SBarry Smith #endif
2556