xref: /petsc/src/ts/interface/ts.c (revision 2bd2b0e6efd31c4ff81d5cb633a913ae833c5993)
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 {
73ace3abfcSBarry Smith   PetscBool      opt,flg;
74dfbe8321SBarry Smith   PetscErrorCode ierr;
75649052a6SBarry Smith   PetscViewer    monviewer;
76eabae89aSBarry Smith   char           monfilename[PETSC_MAX_PATH_LEN];
77089b2837SJed Brown   SNES           snes;
78bdad233fSMatthew Knepley 
79bdad233fSMatthew Knepley   PetscFunctionBegin;
800700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
813194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr);
82a43b19c4SJed Brown     /* Handle TS type options */
83a43b19c4SJed Brown     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
84bdad233fSMatthew Knepley 
85bdad233fSMatthew Knepley     /* Handle generic TS options */
86bdad233fSMatthew Knepley     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
87bdad233fSMatthew Knepley     ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
88d8cd7023SBarry Smith     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);CHKERRQ(ierr);
89d8cd7023SBarry Smith     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&ts->time_step,PETSC_NULL);CHKERRQ(ierr);
90a43b19c4SJed Brown     opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time;
91a43b19c4SJed Brown     ierr = PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);CHKERRQ(ierr);
92461e00b3SSean Farley     if (flg) {ierr = TSSetExactFinalTime(ts,opt);CHKERRQ(ierr);}
93193ac0bcSJed 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);
94193ac0bcSJed Brown     ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections","",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr);
95193ac0bcSJed Brown     ierr = PetscOptionsBool("-ts_error_if_step_failed","Error if no step succeeds","",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr);
96bdad233fSMatthew Knepley 
97bdad233fSMatthew Knepley     /* Monitor options */
98a6570f20SBarry Smith     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
99eabae89aSBarry Smith     if (flg) {
100649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr);
101649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
102bdad233fSMatthew Knepley     }
1035180491cSLisandro Dalcin     ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1045180491cSLisandro Dalcin     if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);}
1055180491cSLisandro Dalcin 
10690d69ab7SBarry Smith     opt  = PETSC_FALSE;
107acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
108a7cc72afSBarry Smith     if (opt) {
109a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
110bdad233fSMatthew Knepley     }
11190d69ab7SBarry Smith     opt  = PETSC_FALSE;
112acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
113a7cc72afSBarry Smith     if (opt) {
1146083293cSBarry Smith       void *ctx;
1156083293cSBarry Smith       ierr = TSMonitorSolutionCreate(ts,PETSC_NULL,&ctx);CHKERRQ(ierr);
1166083293cSBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);CHKERRQ(ierr);
117bdad233fSMatthew Knepley     }
118bdad233fSMatthew Knepley 
119d52bd9f3SBarry Smith     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
120d52bd9f3SBarry Smith     if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
121d52bd9f3SBarry Smith 
122bdad233fSMatthew Knepley     /* Handle specific TS options */
123abc0a331SBarry Smith     if (ts->ops->setfromoptions) {
124bdad233fSMatthew Knepley       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
125bdad233fSMatthew Knepley     }
1265d973c19SBarry Smith 
1275d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
1285d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
129bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
130bdad233fSMatthew Knepley   PetscFunctionReturn(0);
131bdad233fSMatthew Knepley }
132bdad233fSMatthew Knepley 
133bdad233fSMatthew Knepley #undef __FUNCT__
134cdcf91faSSean Farley #undef __FUNCT__
1354a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
136a7a1495cSBarry Smith /*@
1378c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
138a7a1495cSBarry Smith       set with TSSetRHSJacobian().
139a7a1495cSBarry Smith 
140a7a1495cSBarry Smith    Collective on TS and Vec
141a7a1495cSBarry Smith 
142a7a1495cSBarry Smith    Input Parameters:
143316643e7SJed Brown +  ts - the TS context
144a7a1495cSBarry Smith .  t - current timestep
145a7a1495cSBarry Smith -  x - input vector
146a7a1495cSBarry Smith 
147a7a1495cSBarry Smith    Output Parameters:
148a7a1495cSBarry Smith +  A - Jacobian matrix
149a7a1495cSBarry Smith .  B - optional preconditioning matrix
150a7a1495cSBarry Smith -  flag - flag indicating matrix structure
151a7a1495cSBarry Smith 
152a7a1495cSBarry Smith    Notes:
153a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
154a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
155a7a1495cSBarry Smith 
15694b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
157a7a1495cSBarry Smith    flag parameter.
158a7a1495cSBarry Smith 
159a7a1495cSBarry Smith    Level: developer
160a7a1495cSBarry Smith 
161a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
162a7a1495cSBarry Smith 
16394b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
164a7a1495cSBarry Smith @*/
1657087cfbeSBarry Smith PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
166a7a1495cSBarry Smith {
167dfbe8321SBarry Smith   PetscErrorCode ierr;
1680e4ef248SJed Brown   PetscInt Xstate;
169a7a1495cSBarry Smith 
170a7a1495cSBarry Smith   PetscFunctionBegin;
1710700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1720700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
173c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
1740e4ef248SJed Brown   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
1750e4ef248SJed Brown   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) {
1760e4ef248SJed Brown     *flg = ts->rhsjacobian.mstructure;
1770e4ef248SJed Brown     PetscFunctionReturn(0);
178f8ede8e7SJed Brown   }
179d90be118SSean Farley 
1803b5f76d0SSean Farley   if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
181d90be118SSean Farley 
1823b5f76d0SSean Farley   if (ts->userops->rhsjacobian) {
183d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
184a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
185a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
1863b5f76d0SSean Farley     ierr = (*ts->userops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
187a7a1495cSBarry Smith     PetscStackPop;
188d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
189a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
1900700a824SBarry Smith     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
1910700a824SBarry Smith     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
192ef66eb69SBarry Smith   } else {
193214bc6a2SJed Brown     ierr = MatZeroEntries(*A);CHKERRQ(ierr);
194214bc6a2SJed Brown     if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);}
195214bc6a2SJed Brown     *flg = SAME_NONZERO_PATTERN;
196ef66eb69SBarry Smith   }
1970e4ef248SJed Brown   ts->rhsjacobian.time = t;
1980e4ef248SJed Brown   ts->rhsjacobian.X = X;
1990e4ef248SJed Brown   ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr);
2000e4ef248SJed Brown   ts->rhsjacobian.mstructure = *flg;
201a7a1495cSBarry Smith   PetscFunctionReturn(0);
202a7a1495cSBarry Smith }
203a7a1495cSBarry Smith 
2044a2ae208SSatish Balay #undef __FUNCT__
2054a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
206316643e7SJed Brown /*@
207d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
208d763cef2SBarry Smith 
209316643e7SJed Brown    Collective on TS and Vec
210316643e7SJed Brown 
211316643e7SJed Brown    Input Parameters:
212316643e7SJed Brown +  ts - the TS context
213316643e7SJed Brown .  t - current time
214316643e7SJed Brown -  x - state vector
215316643e7SJed Brown 
216316643e7SJed Brown    Output Parameter:
217316643e7SJed Brown .  y - right hand side
218316643e7SJed Brown 
219316643e7SJed Brown    Note:
220316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
221316643e7SJed Brown    is used internally within the nonlinear solvers.
222316643e7SJed Brown 
223316643e7SJed Brown    Level: developer
224316643e7SJed Brown 
225316643e7SJed Brown .keywords: TS, compute
226316643e7SJed Brown 
227316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction()
228316643e7SJed Brown @*/
229dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
230d763cef2SBarry Smith {
231dfbe8321SBarry Smith   PetscErrorCode ierr;
232d763cef2SBarry Smith 
233d763cef2SBarry Smith   PetscFunctionBegin;
2340700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2350700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2360700a824SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
237d763cef2SBarry Smith 
2383b5f76d0SSean Farley   if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
239d763cef2SBarry Smith 
240d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
2413b5f76d0SSean Farley   if (ts->userops->rhsfunction) {
242d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
2433b5f76d0SSean Farley     ierr = (*ts->userops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
244d763cef2SBarry Smith     PetscStackPop;
245214bc6a2SJed Brown   } else {
246214bc6a2SJed Brown     ierr = VecZeroEntries(y);CHKERRQ(ierr);
247b2cd27e8SJed Brown   }
24844a41b28SSean Farley 
249d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
250d763cef2SBarry Smith   PetscFunctionReturn(0);
251d763cef2SBarry Smith }
252d763cef2SBarry Smith 
2534a2ae208SSatish Balay #undef __FUNCT__
254214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSVec_Private"
255214bc6a2SJed Brown static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
256214bc6a2SJed Brown {
2572dd45cf8SJed Brown   Vec            F;
258214bc6a2SJed Brown   PetscErrorCode ierr;
259214bc6a2SJed Brown 
260214bc6a2SJed Brown   PetscFunctionBegin;
2610da9a4f0SJed Brown   *Frhs = PETSC_NULL;
2622dd45cf8SJed Brown   ierr = TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
263214bc6a2SJed Brown   if (!ts->Frhs) {
2642dd45cf8SJed Brown     ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr);
265214bc6a2SJed Brown   }
266214bc6a2SJed Brown   *Frhs = ts->Frhs;
267214bc6a2SJed Brown   PetscFunctionReturn(0);
268214bc6a2SJed Brown }
269214bc6a2SJed Brown 
270214bc6a2SJed Brown #undef __FUNCT__
271214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSMats_Private"
272214bc6a2SJed Brown static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
273214bc6a2SJed Brown {
274214bc6a2SJed Brown   Mat            A,B;
2752dd45cf8SJed Brown   PetscErrorCode ierr;
276214bc6a2SJed Brown 
277214bc6a2SJed Brown   PetscFunctionBegin;
278214bc6a2SJed Brown   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
279214bc6a2SJed Brown   if (Arhs) {
280214bc6a2SJed Brown     if (!ts->Arhs) {
281214bc6a2SJed Brown       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
282214bc6a2SJed Brown     }
283214bc6a2SJed Brown     *Arhs = ts->Arhs;
284214bc6a2SJed Brown   }
285214bc6a2SJed Brown   if (Brhs) {
286214bc6a2SJed Brown     if (!ts->Brhs) {
287214bc6a2SJed Brown       ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
288214bc6a2SJed Brown     }
289214bc6a2SJed Brown     *Brhs = ts->Brhs;
290214bc6a2SJed Brown   }
291214bc6a2SJed Brown   PetscFunctionReturn(0);
292214bc6a2SJed Brown }
293214bc6a2SJed Brown 
294214bc6a2SJed Brown #undef __FUNCT__
295316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction"
296316643e7SJed Brown /*@
297316643e7SJed Brown    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
298316643e7SJed Brown 
299316643e7SJed Brown    Collective on TS and Vec
300316643e7SJed Brown 
301316643e7SJed Brown    Input Parameters:
302316643e7SJed Brown +  ts - the TS context
303316643e7SJed Brown .  t - current time
304316643e7SJed Brown .  X - state vector
305214bc6a2SJed Brown .  Xdot - time derivative of state vector
306214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
307316643e7SJed Brown 
308316643e7SJed Brown    Output Parameter:
309316643e7SJed Brown .  Y - right hand side
310316643e7SJed Brown 
311316643e7SJed Brown    Note:
312316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
313316643e7SJed Brown    is used internally within the nonlinear solvers.
314316643e7SJed Brown 
315316643e7SJed Brown    If the user did did not write their equations in implicit form, this
316316643e7SJed Brown    function recasts them in implicit form.
317316643e7SJed Brown 
318316643e7SJed Brown    Level: developer
319316643e7SJed Brown 
320316643e7SJed Brown .keywords: TS, compute
321316643e7SJed Brown 
322316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction()
323316643e7SJed Brown @*/
324214bc6a2SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
325316643e7SJed Brown {
326316643e7SJed Brown   PetscErrorCode ierr;
327316643e7SJed Brown 
328316643e7SJed Brown   PetscFunctionBegin;
3290700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3300700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
3310700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
3320700a824SBarry Smith   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
333316643e7SJed Brown 
334d90be118SSean Farley   if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
335d90be118SSean Farley 
336316643e7SJed Brown   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
3373b5f76d0SSean Farley   if (ts->userops->ifunction) {
338316643e7SJed Brown     PetscStackPush("TS user implicit function");
3393b5f76d0SSean Farley     ierr = (*ts->userops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr);
340316643e7SJed Brown     PetscStackPop;
341214bc6a2SJed Brown   }
342214bc6a2SJed Brown   if (imex) {
3433b5f76d0SSean Farley     if (!ts->userops->ifunction) {
3442dd45cf8SJed Brown       ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);
3452dd45cf8SJed Brown     }
3462dd45cf8SJed Brown   } else if (ts->userops->rhsfunction) {
3472dd45cf8SJed Brown     if (ts->userops->ifunction) {
348214bc6a2SJed Brown       Vec Frhs;
349214bc6a2SJed Brown       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
350214bc6a2SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr);
351214bc6a2SJed Brown       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
3522dd45cf8SJed Brown     } else {
3532dd45cf8SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr);
3542dd45cf8SJed Brown       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
355316643e7SJed Brown     }
3564a6899ffSJed Brown   }
357316643e7SJed Brown   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
358316643e7SJed Brown   PetscFunctionReturn(0);
359316643e7SJed Brown }
360316643e7SJed Brown 
361316643e7SJed Brown #undef __FUNCT__
362316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian"
363316643e7SJed Brown /*@
364316643e7SJed Brown    TSComputeIJacobian - Evaluates the Jacobian of the DAE
365316643e7SJed Brown 
366316643e7SJed Brown    Collective on TS and Vec
367316643e7SJed Brown 
368316643e7SJed Brown    Input
369316643e7SJed Brown       Input Parameters:
370316643e7SJed Brown +  ts - the TS context
371316643e7SJed Brown .  t - current timestep
372316643e7SJed Brown .  X - state vector
373316643e7SJed Brown .  Xdot - time derivative of state vector
374214bc6a2SJed Brown .  shift - shift to apply, see note below
375214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
376316643e7SJed Brown 
377316643e7SJed Brown    Output Parameters:
378316643e7SJed Brown +  A - Jacobian matrix
379316643e7SJed Brown .  B - optional preconditioning matrix
380316643e7SJed Brown -  flag - flag indicating matrix structure
381316643e7SJed Brown 
382316643e7SJed Brown    Notes:
383316643e7SJed Brown    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
384316643e7SJed Brown 
385316643e7SJed Brown    dF/dX + shift*dF/dXdot
386316643e7SJed Brown 
387316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
388316643e7SJed Brown    is used internally within the nonlinear solvers.
389316643e7SJed Brown 
390316643e7SJed Brown    Level: developer
391316643e7SJed Brown 
392316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix
393316643e7SJed Brown 
394316643e7SJed Brown .seealso:  TSSetIJacobian()
39563495f91SJed Brown @*/
396214bc6a2SJed Brown PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
397316643e7SJed Brown {
3980026cea9SSean Farley   PetscInt Xstate, Xdotstate;
399316643e7SJed Brown   PetscErrorCode ierr;
400316643e7SJed Brown 
401316643e7SJed Brown   PetscFunctionBegin;
4020700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4030700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
4040700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
405316643e7SJed Brown   PetscValidPointer(A,6);
4060700a824SBarry Smith   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
407316643e7SJed Brown   PetscValidPointer(B,7);
4080700a824SBarry Smith   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
409316643e7SJed Brown   PetscValidPointer(flg,8);
4100026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
4110026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr);
4120026cea9SSean Farley   if (ts->ijacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->ijacobian.X == X && ts->ijacobian.Xstate == Xstate && ts->ijacobian.Xdot == Xdot && ts->ijacobian.Xdotstate == Xdotstate && ts->ijacobian.imex == imex))) {
4130026cea9SSean Farley     *flg = ts->ijacobian.mstructure;
4140026cea9SSean Farley     ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr);
4150026cea9SSean Farley     PetscFunctionReturn(0);
4160026cea9SSean Farley   }
417316643e7SJed Brown 
4183b5f76d0SSean Farley   if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
419316643e7SJed Brown 
4204e684422SJed Brown   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
421316643e7SJed Brown   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
4223b5f76d0SSean Farley   if (ts->userops->ijacobian) {
4232dd45cf8SJed Brown     *flg = DIFFERENT_NONZERO_PATTERN;
424316643e7SJed Brown     PetscStackPush("TS user implicit Jacobian");
4253b5f76d0SSean Farley     ierr = (*ts->userops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr);
426316643e7SJed Brown     PetscStackPop;
427214bc6a2SJed Brown     /* make sure user returned a correct Jacobian and preconditioner */
428214bc6a2SJed Brown     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
429214bc6a2SJed Brown     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
4304a6899ffSJed Brown   }
431214bc6a2SJed Brown   if (imex) {
4323b5f76d0SSean Farley     if (!ts->userops->ijacobian) {  /* system was written as Xdot = F(t,X) */
433214bc6a2SJed Brown       ierr = MatZeroEntries(*A);CHKERRQ(ierr);
4342dd45cf8SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
435214bc6a2SJed Brown       if (*A != *B) {
436214bc6a2SJed Brown         ierr = MatZeroEntries(*B);CHKERRQ(ierr);
437214bc6a2SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
438214bc6a2SJed Brown       }
439214bc6a2SJed Brown       *flg = SAME_PRECONDITIONER;
440214bc6a2SJed Brown     }
441214bc6a2SJed Brown   } else {
4423b5f76d0SSean Farley     if (!ts->userops->ijacobian) {
443214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr);
444214bc6a2SJed Brown       ierr = MatScale(*A,-1);CHKERRQ(ierr);
445214bc6a2SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
446316643e7SJed Brown       if (*A != *B) {
447316643e7SJed Brown         ierr = MatScale(*B,-1);CHKERRQ(ierr);
448316643e7SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
449316643e7SJed Brown       }
4503b5f76d0SSean Farley     } else if (ts->userops->rhsjacobian) {
451214bc6a2SJed Brown       Mat Arhs,Brhs;
452214bc6a2SJed Brown       MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
453214bc6a2SJed Brown       ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
454214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
455214bc6a2SJed Brown       axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
456214bc6a2SJed Brown       ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr);
457214bc6a2SJed Brown       if (*A != *B) {
458214bc6a2SJed Brown         ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr);
459214bc6a2SJed Brown       }
460214bc6a2SJed Brown       *flg = PetscMin(*flg,flg2);
461214bc6a2SJed Brown     }
462316643e7SJed Brown   }
4630026cea9SSean Farley 
4640026cea9SSean Farley   ts->ijacobian.time = t;
4650026cea9SSean Farley   ts->ijacobian.X = X;
4660026cea9SSean Farley   ts->ijacobian.Xdot = Xdot;
4670026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr);
4680026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr);
4690026cea9SSean Farley   ts->ijacobian.shift = shift;
4700026cea9SSean Farley   ts->ijacobian.imex = imex;
4710026cea9SSean Farley   ts->ijacobian.mstructure = *flg;
472316643e7SJed Brown   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
473316643e7SJed Brown   PetscFunctionReturn(0);
474316643e7SJed Brown }
475316643e7SJed Brown 
476316643e7SJed Brown #undef __FUNCT__
4774a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
478d763cef2SBarry Smith /*@C
479d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
480d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
481d763cef2SBarry Smith 
4823f9fe445SBarry Smith     Logically Collective on TS
483d763cef2SBarry Smith 
484d763cef2SBarry Smith     Input Parameters:
485d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
486ca94891dSJed Brown .   r - vector to put the computed right hand side (or PETSC_NULL to have it created)
487d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
488d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
489d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
490d763cef2SBarry Smith 
491d763cef2SBarry Smith     Calling sequence of func:
49287828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
493d763cef2SBarry Smith 
494d763cef2SBarry Smith +   t - current timestep
495d763cef2SBarry Smith .   u - input vector
496d763cef2SBarry Smith .   F - function vector
497d763cef2SBarry Smith -   ctx - [optional] user-defined function context
498d763cef2SBarry Smith 
499d763cef2SBarry Smith     Important:
50076f2fa84SHong Zhang     The user MUST call either this routine or TSSetMatrices().
501d763cef2SBarry Smith 
502d763cef2SBarry Smith     Level: beginner
503d763cef2SBarry Smith 
504d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
505d763cef2SBarry Smith 
50676f2fa84SHong Zhang .seealso: TSSetMatrices()
507d763cef2SBarry Smith @*/
508089b2837SJed Brown PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
509d763cef2SBarry Smith {
510089b2837SJed Brown   PetscErrorCode ierr;
511089b2837SJed Brown   SNES           snes;
512d763cef2SBarry Smith 
513089b2837SJed Brown   PetscFunctionBegin;
5140700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
515ca94891dSJed Brown   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
5163b5f76d0SSean Farley   if (f)   ts->userops->rhsfunction = f;
5174e684422SJed Brown   if (ctx) ts->funP                 = ctx;
518089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
519089b2837SJed Brown   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
520d763cef2SBarry Smith   PetscFunctionReturn(0);
521d763cef2SBarry Smith }
522d763cef2SBarry Smith 
5234a2ae208SSatish Balay #undef __FUNCT__
5244a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
525d763cef2SBarry Smith /*@C
526d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
527d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
52876f2fa84SHong Zhang    Use TSSetMatrices() for linear problems.
529d763cef2SBarry Smith 
5303f9fe445SBarry Smith    Logically Collective on TS
531d763cef2SBarry Smith 
532d763cef2SBarry Smith    Input Parameters:
533d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
534d763cef2SBarry Smith .  A   - Jacobian matrix
535d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
536d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
537d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
538d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
539d763cef2SBarry Smith 
540d763cef2SBarry Smith    Calling sequence of func:
54187828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
542d763cef2SBarry Smith 
543d763cef2SBarry Smith +  t - current timestep
544d763cef2SBarry Smith .  u - input vector
545d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
546d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
547d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
54894b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
549d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
550d763cef2SBarry Smith 
551d763cef2SBarry Smith    Notes:
55294b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
553d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
554d763cef2SBarry Smith 
555d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
556d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
55756335db2SHong Zhang    completely new matrix structure (not just different matrix elements)
558d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
559d763cef2SBarry Smith    throughout the global iterations.
560d763cef2SBarry Smith 
561d763cef2SBarry Smith    Level: beginner
562d763cef2SBarry Smith 
563d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
564d763cef2SBarry Smith 
565d763cef2SBarry Smith .seealso: TSDefaultComputeJacobianColor(),
56676f2fa84SHong Zhang           SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices()
567d763cef2SBarry Smith 
568d763cef2SBarry Smith @*/
569089b2837SJed Brown PetscErrorCode  TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
570d763cef2SBarry Smith {
571277b19d0SLisandro Dalcin   PetscErrorCode ierr;
572089b2837SJed Brown   SNES           snes;
573277b19d0SLisandro Dalcin 
574d763cef2SBarry Smith   PetscFunctionBegin;
5750700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
576277b19d0SLisandro Dalcin   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
577277b19d0SLisandro Dalcin   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
578277b19d0SLisandro Dalcin   if (A) PetscCheckSameComm(ts,1,A,2);
579277b19d0SLisandro Dalcin   if (B) PetscCheckSameComm(ts,1,B,3);
580d763cef2SBarry Smith 
5813b5f76d0SSean Farley   if (f)   ts->userops->rhsjacobian = f;
582277b19d0SLisandro Dalcin   if (ctx) ts->jacP                 = ctx;
583089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
584dabe61deSJed Brown   if (!ts->userops->ijacobian) {
585089b2837SJed Brown     ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
5860e4ef248SJed Brown   }
5870e4ef248SJed Brown   if (A) {
5880e4ef248SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
5890e4ef248SJed Brown     ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
5900e4ef248SJed Brown     ts->Arhs = A;
5910e4ef248SJed Brown   }
5920e4ef248SJed Brown   if (B) {
5930e4ef248SJed Brown     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
5940e4ef248SJed Brown     ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
5950e4ef248SJed Brown     ts->Brhs = B;
5960e4ef248SJed Brown   }
597d763cef2SBarry Smith   PetscFunctionReturn(0);
598d763cef2SBarry Smith }
599d763cef2SBarry Smith 
600316643e7SJed Brown 
601316643e7SJed Brown #undef __FUNCT__
602316643e7SJed Brown #define __FUNCT__ "TSSetIFunction"
603316643e7SJed Brown /*@C
604316643e7SJed Brown    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
605316643e7SJed Brown 
6063f9fe445SBarry Smith    Logically Collective on TS
607316643e7SJed Brown 
608316643e7SJed Brown    Input Parameters:
609316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
610ca94891dSJed Brown .  r   - vector to hold the residual (or PETSC_NULL to have it created internally)
611316643e7SJed Brown .  f   - the function evaluation routine
612316643e7SJed Brown -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
613316643e7SJed Brown 
614316643e7SJed Brown    Calling sequence of f:
615316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
616316643e7SJed Brown 
617316643e7SJed Brown +  t   - time at step/stage being solved
618316643e7SJed Brown .  u   - state vector
619316643e7SJed Brown .  u_t - time derivative of state vector
620316643e7SJed Brown .  F   - function vector
621316643e7SJed Brown -  ctx - [optional] user-defined context for matrix evaluation routine
622316643e7SJed Brown 
623316643e7SJed Brown    Important:
624316643e7SJed Brown    The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices().  This routine must be used when not solving an ODE.
625316643e7SJed Brown 
626316643e7SJed Brown    Level: beginner
627316643e7SJed Brown 
628316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian
629316643e7SJed Brown 
630316643e7SJed Brown .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian()
631316643e7SJed Brown @*/
632089b2837SJed Brown PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
633316643e7SJed Brown {
634089b2837SJed Brown   PetscErrorCode ierr;
635089b2837SJed Brown   SNES           snes;
636316643e7SJed Brown 
637316643e7SJed Brown   PetscFunctionBegin;
6380700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
639ca94891dSJed Brown   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
6403b5f76d0SSean Farley   if (f)   ts->userops->ifunction = f;
641089b2837SJed Brown   if (ctx) ts->funP           = ctx;
642089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
643089b2837SJed Brown   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
644089b2837SJed Brown   PetscFunctionReturn(0);
645089b2837SJed Brown }
646089b2837SJed Brown 
647089b2837SJed Brown #undef __FUNCT__
648089b2837SJed Brown #define __FUNCT__ "TSGetIFunction"
649089b2837SJed Brown /*@C
650089b2837SJed Brown    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
651089b2837SJed Brown 
652089b2837SJed Brown    Not Collective
653089b2837SJed Brown 
654089b2837SJed Brown    Input Parameter:
655089b2837SJed Brown .  ts - the TS context
656089b2837SJed Brown 
657089b2837SJed Brown    Output Parameter:
658089b2837SJed Brown +  r - vector to hold residual (or PETSC_NULL)
659089b2837SJed Brown .  func - the function to compute residual (or PETSC_NULL)
660089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
661089b2837SJed Brown 
662089b2837SJed Brown    Level: advanced
663089b2837SJed Brown 
664089b2837SJed Brown .keywords: TS, nonlinear, get, function
665089b2837SJed Brown 
666089b2837SJed Brown .seealso: TSSetIFunction(), SNESGetFunction()
667089b2837SJed Brown @*/
668089b2837SJed Brown PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
669089b2837SJed Brown {
670089b2837SJed Brown   PetscErrorCode ierr;
671089b2837SJed Brown   SNES snes;
672089b2837SJed Brown 
673089b2837SJed Brown   PetscFunctionBegin;
674089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
675089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
676089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
6773b5f76d0SSean Farley   if (func) *func = ts->userops->ifunction;
678089b2837SJed Brown   if (ctx)  *ctx  = ts->funP;
679089b2837SJed Brown   PetscFunctionReturn(0);
680089b2837SJed Brown }
681089b2837SJed Brown 
682089b2837SJed Brown #undef __FUNCT__
683089b2837SJed Brown #define __FUNCT__ "TSGetRHSFunction"
684089b2837SJed Brown /*@C
685089b2837SJed Brown    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
686089b2837SJed Brown 
687089b2837SJed Brown    Not Collective
688089b2837SJed Brown 
689089b2837SJed Brown    Input Parameter:
690089b2837SJed Brown .  ts - the TS context
691089b2837SJed Brown 
692089b2837SJed Brown    Output Parameter:
693089b2837SJed Brown +  r - vector to hold computed right hand side (or PETSC_NULL)
694089b2837SJed Brown .  func - the function to compute right hand side (or PETSC_NULL)
695089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
696089b2837SJed Brown 
697089b2837SJed Brown    Level: advanced
698089b2837SJed Brown 
699089b2837SJed Brown .keywords: TS, nonlinear, get, function
700089b2837SJed Brown 
701089b2837SJed Brown .seealso: TSSetRhsfunction(), SNESGetFunction()
702089b2837SJed Brown @*/
703089b2837SJed Brown PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
704089b2837SJed Brown {
705089b2837SJed Brown   PetscErrorCode ierr;
706089b2837SJed Brown   SNES snes;
707089b2837SJed Brown 
708089b2837SJed Brown   PetscFunctionBegin;
709089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
710089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
711089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
7123b5f76d0SSean Farley   if (func) *func = ts->userops->rhsfunction;
713089b2837SJed Brown   if (ctx)  *ctx  = ts->funP;
714316643e7SJed Brown   PetscFunctionReturn(0);
715316643e7SJed Brown }
716316643e7SJed Brown 
717316643e7SJed Brown #undef __FUNCT__
718316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian"
719316643e7SJed Brown /*@C
720a4f0a591SBarry Smith    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
721a4f0a591SBarry Smith         you provided with TSSetIFunction().
722316643e7SJed Brown 
7233f9fe445SBarry Smith    Logically Collective on TS
724316643e7SJed Brown 
725316643e7SJed Brown    Input Parameters:
726316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
727316643e7SJed Brown .  A   - Jacobian matrix
728316643e7SJed Brown .  B   - preconditioning matrix for A (may be same as A)
729316643e7SJed Brown .  f   - the Jacobian evaluation routine
730316643e7SJed Brown -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
731316643e7SJed Brown 
732316643e7SJed Brown    Calling sequence of f:
7331b4a444bSJed Brown $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
734316643e7SJed Brown 
735316643e7SJed Brown +  t    - time at step/stage being solved
7361b4a444bSJed Brown .  U    - state vector
7371b4a444bSJed Brown .  U_t  - time derivative of state vector
738316643e7SJed Brown .  a    - shift
7391b4a444bSJed Brown .  A    - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
740316643e7SJed Brown .  B    - preconditioning matrix for A, may be same as A
741316643e7SJed Brown .  flag - flag indicating information about the preconditioner matrix
742316643e7SJed Brown           structure (same as flag in KSPSetOperators())
743316643e7SJed Brown -  ctx  - [optional] user-defined context for matrix evaluation routine
744316643e7SJed Brown 
745316643e7SJed Brown    Notes:
746316643e7SJed Brown    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
747316643e7SJed Brown 
748a4f0a591SBarry Smith    The matrix dF/dU + a*dF/dU_t you provide turns out to be
749a4f0a591SBarry 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.
750a4f0a591SBarry Smith    The time integrator internally approximates U_t by W+a*U where the positive "shift"
751a4f0a591SBarry Smith    a and vector W depend on the integration method, step size, and past states. For example with
752a4f0a591SBarry Smith    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
753a4f0a591SBarry Smith    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
754a4f0a591SBarry Smith 
755316643e7SJed Brown    Level: beginner
756316643e7SJed Brown 
757316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian
758316643e7SJed Brown 
759316643e7SJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian()
760316643e7SJed Brown 
761316643e7SJed Brown @*/
7627087cfbeSBarry Smith PetscErrorCode  TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
763316643e7SJed Brown {
764316643e7SJed Brown   PetscErrorCode ierr;
765089b2837SJed Brown   SNES           snes;
766316643e7SJed Brown 
767316643e7SJed Brown   PetscFunctionBegin;
7680700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7690700a824SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
7700700a824SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
771316643e7SJed Brown   if (A) PetscCheckSameComm(ts,1,A,2);
772316643e7SJed Brown   if (B) PetscCheckSameComm(ts,1,B,3);
7733b5f76d0SSean Farley   if (f)   ts->userops->ijacobian = f;
774316643e7SJed Brown   if (ctx) ts->jacP           = ctx;
775089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
776089b2837SJed Brown   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
777316643e7SJed Brown   PetscFunctionReturn(0);
778316643e7SJed Brown }
779316643e7SJed Brown 
7804a2ae208SSatish Balay #undef __FUNCT__
7814a2ae208SSatish Balay #define __FUNCT__ "TSView"
7827e2c5f70SBarry Smith /*@C
783d763cef2SBarry Smith     TSView - Prints the TS data structure.
784d763cef2SBarry Smith 
7854c49b128SBarry Smith     Collective on TS
786d763cef2SBarry Smith 
787d763cef2SBarry Smith     Input Parameters:
788d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
789d763cef2SBarry Smith -   viewer - visualization context
790d763cef2SBarry Smith 
791d763cef2SBarry Smith     Options Database Key:
792d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
793d763cef2SBarry Smith 
794d763cef2SBarry Smith     Notes:
795d763cef2SBarry Smith     The available visualization contexts include
796b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
797b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
798d763cef2SBarry Smith          output where only the first processor opens
799d763cef2SBarry Smith          the file.  All other processors send their
800d763cef2SBarry Smith          data to the first processor to print.
801d763cef2SBarry Smith 
802d763cef2SBarry Smith     The user can open an alternative visualization context with
803b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
804d763cef2SBarry Smith 
805d763cef2SBarry Smith     Level: beginner
806d763cef2SBarry Smith 
807d763cef2SBarry Smith .keywords: TS, timestep, view
808d763cef2SBarry Smith 
809b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
810d763cef2SBarry Smith @*/
8117087cfbeSBarry Smith PetscErrorCode  TSView(TS ts,PetscViewer viewer)
812d763cef2SBarry Smith {
813dfbe8321SBarry Smith   PetscErrorCode ierr;
814a313700dSBarry Smith   const TSType   type;
815089b2837SJed Brown   PetscBool      iascii,isstring,isundials;
816d763cef2SBarry Smith 
817d763cef2SBarry Smith   PetscFunctionBegin;
8180700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8193050cee2SBarry Smith   if (!viewer) {
8207adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
8213050cee2SBarry Smith   }
8220700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
823c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
824fd16b177SBarry Smith 
8252692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8262692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
82732077d6dSBarry Smith   if (iascii) {
828317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr);
82977431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
830a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
831d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
83277431f27SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr);
833c610991cSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr);
834d763cef2SBarry Smith     }
83577431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr);
836193ac0bcSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
837d52bd9f3SBarry Smith     if (ts->ops->view) {
838d52bd9f3SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
839d52bd9f3SBarry Smith       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
840d52bd9f3SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
841d52bd9f3SBarry Smith     }
8420f5bd95cSBarry Smith   } else if (isstring) {
843a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
844b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
845d763cef2SBarry Smith   }
846b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
847089b2837SJed Brown   ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
848b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
849d763cef2SBarry Smith   PetscFunctionReturn(0);
850d763cef2SBarry Smith }
851d763cef2SBarry Smith 
852d763cef2SBarry Smith 
8534a2ae208SSatish Balay #undef __FUNCT__
8544a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
855b07ff414SBarry Smith /*@
856d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
857d763cef2SBarry Smith    the timesteppers.
858d763cef2SBarry Smith 
8593f9fe445SBarry Smith    Logically Collective on TS
860d763cef2SBarry Smith 
861d763cef2SBarry Smith    Input Parameters:
862d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
863d763cef2SBarry Smith -  usrP - optional user context
864d763cef2SBarry Smith 
865d763cef2SBarry Smith    Level: intermediate
866d763cef2SBarry Smith 
867d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
868d763cef2SBarry Smith 
869d763cef2SBarry Smith .seealso: TSGetApplicationContext()
870d763cef2SBarry Smith @*/
8717087cfbeSBarry Smith PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
872d763cef2SBarry Smith {
873d763cef2SBarry Smith   PetscFunctionBegin;
8740700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
875d763cef2SBarry Smith   ts->user = usrP;
876d763cef2SBarry Smith   PetscFunctionReturn(0);
877d763cef2SBarry Smith }
878d763cef2SBarry Smith 
8794a2ae208SSatish Balay #undef __FUNCT__
8804a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
881b07ff414SBarry Smith /*@
882d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
883d763cef2SBarry Smith     timestepper.
884d763cef2SBarry Smith 
885d763cef2SBarry Smith     Not Collective
886d763cef2SBarry Smith 
887d763cef2SBarry Smith     Input Parameter:
888d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
889d763cef2SBarry Smith 
890d763cef2SBarry Smith     Output Parameter:
891d763cef2SBarry Smith .   usrP - user context
892d763cef2SBarry Smith 
893d763cef2SBarry Smith     Level: intermediate
894d763cef2SBarry Smith 
895d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
896d763cef2SBarry Smith 
897d763cef2SBarry Smith .seealso: TSSetApplicationContext()
898d763cef2SBarry Smith @*/
899e71120c6SJed Brown PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
900d763cef2SBarry Smith {
901d763cef2SBarry Smith   PetscFunctionBegin;
9020700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
903e71120c6SJed Brown   *(void**)usrP = ts->user;
904d763cef2SBarry Smith   PetscFunctionReturn(0);
905d763cef2SBarry Smith }
906d763cef2SBarry Smith 
9074a2ae208SSatish Balay #undef __FUNCT__
9084a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
909d763cef2SBarry Smith /*@
910d763cef2SBarry Smith    TSGetTimeStepNumber - Gets the current number of timesteps.
911d763cef2SBarry Smith 
912d763cef2SBarry Smith    Not Collective
913d763cef2SBarry Smith 
914d763cef2SBarry Smith    Input Parameter:
915d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
916d763cef2SBarry Smith 
917d763cef2SBarry Smith    Output Parameter:
918d763cef2SBarry Smith .  iter - number steps so far
919d763cef2SBarry Smith 
920d763cef2SBarry Smith    Level: intermediate
921d763cef2SBarry Smith 
922d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
923d763cef2SBarry Smith @*/
9247087cfbeSBarry Smith PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
925d763cef2SBarry Smith {
926d763cef2SBarry Smith   PetscFunctionBegin;
9270700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
9284482741eSBarry Smith   PetscValidIntPointer(iter,2);
929d763cef2SBarry Smith   *iter = ts->steps;
930d763cef2SBarry Smith   PetscFunctionReturn(0);
931d763cef2SBarry Smith }
932d763cef2SBarry Smith 
9334a2ae208SSatish Balay #undef __FUNCT__
9344a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
935d763cef2SBarry Smith /*@
936d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
937d763cef2SBarry Smith    as well as the initial time.
938d763cef2SBarry Smith 
9393f9fe445SBarry Smith    Logically Collective on TS
940d763cef2SBarry Smith 
941d763cef2SBarry Smith    Input Parameters:
942d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
943d763cef2SBarry Smith .  initial_time - the initial time
944d763cef2SBarry Smith -  time_step - the size of the timestep
945d763cef2SBarry Smith 
946d763cef2SBarry Smith    Level: intermediate
947d763cef2SBarry Smith 
948d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
949d763cef2SBarry Smith 
950d763cef2SBarry Smith .keywords: TS, set, initial, timestep
951d763cef2SBarry Smith @*/
9527087cfbeSBarry Smith PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
953d763cef2SBarry Smith {
954e144a568SJed Brown   PetscErrorCode ierr;
955e144a568SJed Brown 
956d763cef2SBarry Smith   PetscFunctionBegin;
9570700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
958e144a568SJed Brown   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
959d8cd7023SBarry Smith   ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr);
960d763cef2SBarry Smith   PetscFunctionReturn(0);
961d763cef2SBarry Smith }
962d763cef2SBarry Smith 
9634a2ae208SSatish Balay #undef __FUNCT__
9644a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
965d763cef2SBarry Smith /*@
966d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
967d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
968d763cef2SBarry Smith 
9693f9fe445SBarry Smith    Logically Collective on TS
970d763cef2SBarry Smith 
971d763cef2SBarry Smith    Input Parameters:
972d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
973d763cef2SBarry Smith -  time_step - the size of the timestep
974d763cef2SBarry Smith 
975d763cef2SBarry Smith    Level: intermediate
976d763cef2SBarry Smith 
977d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
978d763cef2SBarry Smith 
979d763cef2SBarry Smith .keywords: TS, set, timestep
980d763cef2SBarry Smith @*/
9817087cfbeSBarry Smith PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
982d763cef2SBarry Smith {
983d763cef2SBarry Smith   PetscFunctionBegin;
9840700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
985c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,time_step,2);
986d763cef2SBarry Smith   ts->time_step = time_step;
987d763cef2SBarry Smith   PetscFunctionReturn(0);
988d763cef2SBarry Smith }
989d763cef2SBarry Smith 
9904a2ae208SSatish Balay #undef __FUNCT__
991a43b19c4SJed Brown #define __FUNCT__ "TSSetExactFinalTime"
992a43b19c4SJed Brown /*@
993a43b19c4SJed Brown    TSSetExactFinalTime - Determines whether to interpolate solution to the
994a43b19c4SJed Brown       exact final time requested by the user or just returns it at the final time
995a43b19c4SJed Brown       it computed.
996a43b19c4SJed Brown 
997a43b19c4SJed Brown   Logically Collective on TS
998a43b19c4SJed Brown 
999a43b19c4SJed Brown    Input Parameter:
1000a43b19c4SJed Brown +   ts - the time-step context
1001a43b19c4SJed Brown -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
1002a43b19c4SJed Brown 
1003a43b19c4SJed Brown    Level: beginner
1004a43b19c4SJed Brown 
1005a43b19c4SJed Brown .seealso: TSSetDuration()
1006a43b19c4SJed Brown @*/
1007a43b19c4SJed Brown PetscErrorCode  TSSetExactFinalTime(TS ts,PetscBool flg)
1008a43b19c4SJed Brown {
1009a43b19c4SJed Brown 
1010a43b19c4SJed Brown   PetscFunctionBegin;
1011a43b19c4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1012d8cd7023SBarry Smith   PetscValidLogicalCollectiveBool(ts,flg,2);
1013a43b19c4SJed Brown   ts->exact_final_time = flg;
1014a43b19c4SJed Brown   PetscFunctionReturn(0);
1015a43b19c4SJed Brown }
1016a43b19c4SJed Brown 
1017a43b19c4SJed Brown #undef __FUNCT__
10184a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
1019d763cef2SBarry Smith /*@
1020d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
1021d763cef2SBarry Smith 
1022d763cef2SBarry Smith    Not Collective
1023d763cef2SBarry Smith 
1024d763cef2SBarry Smith    Input Parameter:
1025d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1026d763cef2SBarry Smith 
1027d763cef2SBarry Smith    Output Parameter:
1028d763cef2SBarry Smith .  dt - the current timestep size
1029d763cef2SBarry Smith 
1030d763cef2SBarry Smith    Level: intermediate
1031d763cef2SBarry Smith 
1032d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1033d763cef2SBarry Smith 
1034d763cef2SBarry Smith .keywords: TS, get, timestep
1035d763cef2SBarry Smith @*/
10367087cfbeSBarry Smith PetscErrorCode  TSGetTimeStep(TS ts,PetscReal* dt)
1037d763cef2SBarry Smith {
1038d763cef2SBarry Smith   PetscFunctionBegin;
10390700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10404482741eSBarry Smith   PetscValidDoublePointer(dt,2);
1041d763cef2SBarry Smith   *dt = ts->time_step;
1042d763cef2SBarry Smith   PetscFunctionReturn(0);
1043d763cef2SBarry Smith }
1044d763cef2SBarry Smith 
10454a2ae208SSatish Balay #undef __FUNCT__
10464a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
1047d8e5e3e6SSatish Balay /*@
1048d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
1049d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
1050d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
1051d763cef2SBarry Smith    the solution at the next timestep has been calculated.
1052d763cef2SBarry Smith 
1053d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
1054d763cef2SBarry Smith 
1055d763cef2SBarry Smith    Input Parameter:
1056d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1057d763cef2SBarry Smith 
1058d763cef2SBarry Smith    Output Parameter:
1059d763cef2SBarry Smith .  v - the vector containing the solution
1060d763cef2SBarry Smith 
1061d763cef2SBarry Smith    Level: intermediate
1062d763cef2SBarry Smith 
1063d763cef2SBarry Smith .seealso: TSGetTimeStep()
1064d763cef2SBarry Smith 
1065d763cef2SBarry Smith .keywords: TS, timestep, get, solution
1066d763cef2SBarry Smith @*/
10677087cfbeSBarry Smith PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1068d763cef2SBarry Smith {
1069d763cef2SBarry Smith   PetscFunctionBegin;
10700700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10714482741eSBarry Smith   PetscValidPointer(v,2);
10728737fe31SLisandro Dalcin   *v = ts->vec_sol;
1073d763cef2SBarry Smith   PetscFunctionReturn(0);
1074d763cef2SBarry Smith }
1075d763cef2SBarry Smith 
1076bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
10774a2ae208SSatish Balay #undef __FUNCT__
1078bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
1079d8e5e3e6SSatish Balay /*@
1080bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
1081d763cef2SBarry Smith 
1082bdad233fSMatthew Knepley   Not collective
1083d763cef2SBarry Smith 
1084bdad233fSMatthew Knepley   Input Parameters:
1085bdad233fSMatthew Knepley + ts   - The TS
1086bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1087d763cef2SBarry Smith .vb
1088d763cef2SBarry Smith          U_t = A U
1089d763cef2SBarry Smith          U_t = A(t) U
1090d763cef2SBarry Smith          U_t = F(t,U)
1091d763cef2SBarry Smith .ve
1092d763cef2SBarry Smith 
1093d763cef2SBarry Smith    Level: beginner
1094d763cef2SBarry Smith 
1095bdad233fSMatthew Knepley .keywords: TS, problem type
1096bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1097d763cef2SBarry Smith @*/
10987087cfbeSBarry Smith PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1099a7cc72afSBarry Smith {
11009e2a6581SJed Brown   PetscErrorCode ierr;
11019e2a6581SJed Brown 
1102d763cef2SBarry Smith   PetscFunctionBegin;
11030700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1104bdad233fSMatthew Knepley   ts->problem_type = type;
11059e2a6581SJed Brown   if (type == TS_LINEAR) {
11069e2a6581SJed Brown     SNES snes;
11079e2a6581SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
11089e2a6581SJed Brown     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
11099e2a6581SJed Brown   }
1110d763cef2SBarry Smith   PetscFunctionReturn(0);
1111d763cef2SBarry Smith }
1112d763cef2SBarry Smith 
1113bdad233fSMatthew Knepley #undef __FUNCT__
1114bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
1115bdad233fSMatthew Knepley /*@C
1116bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
1117bdad233fSMatthew Knepley 
1118bdad233fSMatthew Knepley   Not collective
1119bdad233fSMatthew Knepley 
1120bdad233fSMatthew Knepley   Input Parameter:
1121bdad233fSMatthew Knepley . ts   - The TS
1122bdad233fSMatthew Knepley 
1123bdad233fSMatthew Knepley   Output Parameter:
1124bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1125bdad233fSMatthew Knepley .vb
1126089b2837SJed Brown          M U_t = A U
1127089b2837SJed Brown          M(t) U_t = A(t) U
1128bdad233fSMatthew Knepley          U_t = F(t,U)
1129bdad233fSMatthew Knepley .ve
1130bdad233fSMatthew Knepley 
1131bdad233fSMatthew Knepley    Level: beginner
1132bdad233fSMatthew Knepley 
1133bdad233fSMatthew Knepley .keywords: TS, problem type
1134bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1135bdad233fSMatthew Knepley @*/
11367087cfbeSBarry Smith PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1137a7cc72afSBarry Smith {
1138bdad233fSMatthew Knepley   PetscFunctionBegin;
11390700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
11404482741eSBarry Smith   PetscValidIntPointer(type,2);
1141bdad233fSMatthew Knepley   *type = ts->problem_type;
1142bdad233fSMatthew Knepley   PetscFunctionReturn(0);
1143bdad233fSMatthew Knepley }
1144d763cef2SBarry Smith 
11454a2ae208SSatish Balay #undef __FUNCT__
11464a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
1147d763cef2SBarry Smith /*@
1148d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
1149d763cef2SBarry Smith    of a timestepper.
1150d763cef2SBarry Smith 
1151d763cef2SBarry Smith    Collective on TS
1152d763cef2SBarry Smith 
1153d763cef2SBarry Smith    Input Parameter:
1154d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1155d763cef2SBarry Smith 
1156d763cef2SBarry Smith    Notes:
1157d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
1158d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
1159d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
1160d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
1161d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
1162d763cef2SBarry Smith 
1163d763cef2SBarry Smith    Level: advanced
1164d763cef2SBarry Smith 
1165d763cef2SBarry Smith .keywords: TS, timestep, setup
1166d763cef2SBarry Smith 
1167d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
1168d763cef2SBarry Smith @*/
11697087cfbeSBarry Smith PetscErrorCode  TSSetUp(TS ts)
1170d763cef2SBarry Smith {
1171dfbe8321SBarry Smith   PetscErrorCode ierr;
1172d763cef2SBarry Smith 
1173d763cef2SBarry Smith   PetscFunctionBegin;
11740700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1175277b19d0SLisandro Dalcin   if (ts->setupcalled) PetscFunctionReturn(0);
1176277b19d0SLisandro Dalcin 
11777adad957SLisandro Dalcin   if (!((PetscObject)ts)->type_name) {
11789596e0b4SJed Brown     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1179d763cef2SBarry Smith   }
1180a43b19c4SJed Brown   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1181277b19d0SLisandro Dalcin 
1182277b19d0SLisandro Dalcin   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1183277b19d0SLisandro Dalcin 
1184277b19d0SLisandro Dalcin   if (ts->ops->setup) {
1185000e7ae3SMatthew Knepley     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1186277b19d0SLisandro Dalcin   }
1187277b19d0SLisandro Dalcin 
1188277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_TRUE;
1189277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1190277b19d0SLisandro Dalcin }
1191277b19d0SLisandro Dalcin 
1192277b19d0SLisandro Dalcin #undef __FUNCT__
1193277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset"
1194277b19d0SLisandro Dalcin /*@
1195277b19d0SLisandro Dalcin    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1196277b19d0SLisandro Dalcin 
1197277b19d0SLisandro Dalcin    Collective on TS
1198277b19d0SLisandro Dalcin 
1199277b19d0SLisandro Dalcin    Input Parameter:
1200277b19d0SLisandro Dalcin .  ts - the TS context obtained from TSCreate()
1201277b19d0SLisandro Dalcin 
1202277b19d0SLisandro Dalcin    Level: beginner
1203277b19d0SLisandro Dalcin 
1204277b19d0SLisandro Dalcin .keywords: TS, timestep, reset
1205277b19d0SLisandro Dalcin 
1206277b19d0SLisandro Dalcin .seealso: TSCreate(), TSSetup(), TSDestroy()
1207277b19d0SLisandro Dalcin @*/
1208277b19d0SLisandro Dalcin PetscErrorCode  TSReset(TS ts)
1209277b19d0SLisandro Dalcin {
1210277b19d0SLisandro Dalcin   PetscErrorCode ierr;
1211277b19d0SLisandro Dalcin 
1212277b19d0SLisandro Dalcin   PetscFunctionBegin;
1213277b19d0SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1214277b19d0SLisandro Dalcin   if (ts->ops->reset) {
1215277b19d0SLisandro Dalcin     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1216277b19d0SLisandro Dalcin   }
1217277b19d0SLisandro Dalcin   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
12184e684422SJed Brown   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
12194e684422SJed Brown   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1220214bc6a2SJed Brown   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
12216bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
122238637c2eSJed Brown   ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);
1223277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_FALSE;
1224d763cef2SBarry Smith   PetscFunctionReturn(0);
1225d763cef2SBarry Smith }
1226d763cef2SBarry Smith 
12274a2ae208SSatish Balay #undef __FUNCT__
12284a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
1229d8e5e3e6SSatish Balay /*@
1230d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
1231d763cef2SBarry Smith    with TSCreate().
1232d763cef2SBarry Smith 
1233d763cef2SBarry Smith    Collective on TS
1234d763cef2SBarry Smith 
1235d763cef2SBarry Smith    Input Parameter:
1236d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1237d763cef2SBarry Smith 
1238d763cef2SBarry Smith    Level: beginner
1239d763cef2SBarry Smith 
1240d763cef2SBarry Smith .keywords: TS, timestepper, destroy
1241d763cef2SBarry Smith 
1242d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
1243d763cef2SBarry Smith @*/
12446bf464f9SBarry Smith PetscErrorCode  TSDestroy(TS *ts)
1245d763cef2SBarry Smith {
12466849ba73SBarry Smith   PetscErrorCode ierr;
1247d763cef2SBarry Smith 
1248d763cef2SBarry Smith   PetscFunctionBegin;
12496bf464f9SBarry Smith   if (!*ts) PetscFunctionReturn(0);
12506bf464f9SBarry Smith   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
12516bf464f9SBarry Smith   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1252d763cef2SBarry Smith 
12536bf464f9SBarry Smith   ierr = TSReset((*ts));CHKERRQ(ierr);
1254277b19d0SLisandro Dalcin 
1255be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
12566bf464f9SBarry Smith   ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr);
12576bf464f9SBarry Smith   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
12586d4c513bSLisandro Dalcin 
12596bf464f9SBarry Smith   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
12606bf464f9SBarry Smith   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
12616bf464f9SBarry Smith   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
12626d4c513bSLisandro Dalcin 
12633b5f76d0SSean Farley   ierr = PetscFree((*ts)->userops);
12643b5f76d0SSean Farley 
1265a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1266d763cef2SBarry Smith   PetscFunctionReturn(0);
1267d763cef2SBarry Smith }
1268d763cef2SBarry Smith 
12694a2ae208SSatish Balay #undef __FUNCT__
12704a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
1271d8e5e3e6SSatish Balay /*@
1272d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1273d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
1274d763cef2SBarry Smith 
1275d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
1276d763cef2SBarry Smith 
1277d763cef2SBarry Smith    Input Parameter:
1278d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1279d763cef2SBarry Smith 
1280d763cef2SBarry Smith    Output Parameter:
1281d763cef2SBarry Smith .  snes - the nonlinear solver context
1282d763cef2SBarry Smith 
1283d763cef2SBarry Smith    Notes:
1284d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
1285d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
128694b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
1287d763cef2SBarry Smith 
1288d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
1289d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
1290d763cef2SBarry Smith 
1291d763cef2SBarry Smith    Level: beginner
1292d763cef2SBarry Smith 
1293d763cef2SBarry Smith .keywords: timestep, get, SNES
1294d763cef2SBarry Smith @*/
12957087cfbeSBarry Smith PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1296d763cef2SBarry Smith {
1297d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1298d372ba47SLisandro Dalcin 
1299d763cef2SBarry Smith   PetscFunctionBegin;
13000700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13014482741eSBarry Smith   PetscValidPointer(snes,2);
1302d372ba47SLisandro Dalcin   if (!ts->snes) {
1303d372ba47SLisandro Dalcin     ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
1304d372ba47SLisandro Dalcin     ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr);
1305d372ba47SLisandro Dalcin     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
13069e2a6581SJed Brown     if (ts->problem_type == TS_LINEAR) {
13079e2a6581SJed Brown       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
13089e2a6581SJed Brown     }
1309d372ba47SLisandro Dalcin   }
1310d763cef2SBarry Smith   *snes = ts->snes;
1311d763cef2SBarry Smith   PetscFunctionReturn(0);
1312d763cef2SBarry Smith }
1313d763cef2SBarry Smith 
13144a2ae208SSatish Balay #undef __FUNCT__
131594b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
1316d8e5e3e6SSatish Balay /*@
131794b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
1318d763cef2SBarry Smith    a TS (timestepper) context.
1319d763cef2SBarry Smith 
132094b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1321d763cef2SBarry Smith 
1322d763cef2SBarry Smith    Input Parameter:
1323d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1324d763cef2SBarry Smith 
1325d763cef2SBarry Smith    Output Parameter:
132694b7f48cSBarry Smith .  ksp - the nonlinear solver context
1327d763cef2SBarry Smith 
1328d763cef2SBarry Smith    Notes:
132994b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1330d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1331d763cef2SBarry Smith    KSP and PC contexts as well.
1332d763cef2SBarry Smith 
133394b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
133494b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1335d763cef2SBarry Smith 
1336d763cef2SBarry Smith    Level: beginner
1337d763cef2SBarry Smith 
133894b7f48cSBarry Smith .keywords: timestep, get, KSP
1339d763cef2SBarry Smith @*/
13407087cfbeSBarry Smith PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1341d763cef2SBarry Smith {
1342d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1343089b2837SJed Brown   SNES           snes;
1344d372ba47SLisandro Dalcin 
1345d763cef2SBarry Smith   PetscFunctionBegin;
13460700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13474482741eSBarry Smith   PetscValidPointer(ksp,2);
134817186662SBarry Smith   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1349e32f2f54SBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1350089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1351089b2837SJed Brown   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
1352d763cef2SBarry Smith   PetscFunctionReturn(0);
1353d763cef2SBarry Smith }
1354d763cef2SBarry Smith 
1355d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1356d763cef2SBarry Smith 
13574a2ae208SSatish Balay #undef __FUNCT__
1358adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1359adb62b0dSMatthew Knepley /*@
1360adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1361adb62b0dSMatthew Knepley    maximum time for iteration.
1362adb62b0dSMatthew Knepley 
13633f9fe445SBarry Smith    Not Collective
1364adb62b0dSMatthew Knepley 
1365adb62b0dSMatthew Knepley    Input Parameters:
1366adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1367adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1368adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1369adb62b0dSMatthew Knepley 
1370adb62b0dSMatthew Knepley    Level: intermediate
1371adb62b0dSMatthew Knepley 
1372adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1373adb62b0dSMatthew Knepley @*/
13747087cfbeSBarry Smith PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1375adb62b0dSMatthew Knepley {
1376adb62b0dSMatthew Knepley   PetscFunctionBegin;
13770700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1378abc0a331SBarry Smith   if (maxsteps) {
13794482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1380adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1381adb62b0dSMatthew Knepley   }
1382abc0a331SBarry Smith   if (maxtime) {
13834482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1384adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1385adb62b0dSMatthew Knepley   }
1386adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1387adb62b0dSMatthew Knepley }
1388adb62b0dSMatthew Knepley 
1389adb62b0dSMatthew Knepley #undef __FUNCT__
13904a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1391d763cef2SBarry Smith /*@
1392d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1393d763cef2SBarry Smith    maximum time for iteration.
1394d763cef2SBarry Smith 
13953f9fe445SBarry Smith    Logically Collective on TS
1396d763cef2SBarry Smith 
1397d763cef2SBarry Smith    Input Parameters:
1398d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1399d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1400d763cef2SBarry Smith -  maxtime - final time to iterate to
1401d763cef2SBarry Smith 
1402d763cef2SBarry Smith    Options Database Keys:
1403d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
1404d763cef2SBarry Smith .  -ts_max_time <maxtime> - Sets maxtime
1405d763cef2SBarry Smith 
1406d763cef2SBarry Smith    Notes:
1407d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1408d763cef2SBarry Smith 
1409d763cef2SBarry Smith    Level: intermediate
1410d763cef2SBarry Smith 
1411d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1412a43b19c4SJed Brown 
1413a43b19c4SJed Brown .seealso: TSSetExactFinalTime()
1414d763cef2SBarry Smith @*/
14157087cfbeSBarry Smith PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1416d763cef2SBarry Smith {
1417d763cef2SBarry Smith   PetscFunctionBegin;
14180700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1419c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1420c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,maxtime,2);
142139b7ec4bSSean Farley   if (maxsteps >= 0) ts->max_steps = maxsteps;
142239b7ec4bSSean Farley   if (maxtime != PETSC_DEFAULT) ts->max_time  = maxtime;
1423d763cef2SBarry Smith   PetscFunctionReturn(0);
1424d763cef2SBarry Smith }
1425d763cef2SBarry Smith 
14264a2ae208SSatish Balay #undef __FUNCT__
14274a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1428d763cef2SBarry Smith /*@
1429d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1430d763cef2SBarry Smith    for use by the TS routines.
1431d763cef2SBarry Smith 
14323f9fe445SBarry Smith    Logically Collective on TS and Vec
1433d763cef2SBarry Smith 
1434d763cef2SBarry Smith    Input Parameters:
1435d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1436d763cef2SBarry Smith -  x - the solution vector
1437d763cef2SBarry Smith 
1438d763cef2SBarry Smith    Level: beginner
1439d763cef2SBarry Smith 
1440d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1441d763cef2SBarry Smith @*/
14427087cfbeSBarry Smith PetscErrorCode  TSSetSolution(TS ts,Vec x)
1443d763cef2SBarry Smith {
14448737fe31SLisandro Dalcin   PetscErrorCode ierr;
14458737fe31SLisandro Dalcin 
1446d763cef2SBarry Smith   PetscFunctionBegin;
14470700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
14480700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
14498737fe31SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
14506bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
14518737fe31SLisandro Dalcin   ts->vec_sol = x;
1452d763cef2SBarry Smith   PetscFunctionReturn(0);
1453d763cef2SBarry Smith }
1454d763cef2SBarry Smith 
1455e74ef692SMatthew Knepley #undef __FUNCT__
1456e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1457ac226902SBarry Smith /*@C
1458000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
14593f2090d5SJed Brown   called once at the beginning of each time step.
1460000e7ae3SMatthew Knepley 
14613f9fe445SBarry Smith   Logically Collective on TS
1462000e7ae3SMatthew Knepley 
1463000e7ae3SMatthew Knepley   Input Parameters:
1464000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1465000e7ae3SMatthew Knepley - func - The function
1466000e7ae3SMatthew Knepley 
1467000e7ae3SMatthew Knepley   Calling sequence of func:
1468000e7ae3SMatthew Knepley . func (TS ts);
1469000e7ae3SMatthew Knepley 
1470000e7ae3SMatthew Knepley   Level: intermediate
1471000e7ae3SMatthew Knepley 
1472000e7ae3SMatthew Knepley .keywords: TS, timestep
1473000e7ae3SMatthew Knepley @*/
14747087cfbeSBarry Smith PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1475000e7ae3SMatthew Knepley {
1476000e7ae3SMatthew Knepley   PetscFunctionBegin;
14770700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1478000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1479000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1480000e7ae3SMatthew Knepley }
1481000e7ae3SMatthew Knepley 
1482e74ef692SMatthew Knepley #undef __FUNCT__
14833f2090d5SJed Brown #define __FUNCT__ "TSPreStep"
14843f2090d5SJed Brown /*@C
14853f2090d5SJed Brown   TSPreStep - Runs the user-defined pre-step function.
14863f2090d5SJed Brown 
14873f2090d5SJed Brown   Collective on TS
14883f2090d5SJed Brown 
14893f2090d5SJed Brown   Input Parameters:
14903f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
14913f2090d5SJed Brown 
14923f2090d5SJed Brown   Notes:
14933f2090d5SJed Brown   TSPreStep() is typically used within time stepping implementations,
14943f2090d5SJed Brown   so most users would not generally call this routine themselves.
14953f2090d5SJed Brown 
14963f2090d5SJed Brown   Level: developer
14973f2090d5SJed Brown 
14983f2090d5SJed Brown .keywords: TS, timestep
14993f2090d5SJed Brown @*/
15007087cfbeSBarry Smith PetscErrorCode  TSPreStep(TS ts)
15013f2090d5SJed Brown {
15023f2090d5SJed Brown   PetscErrorCode ierr;
15033f2090d5SJed Brown 
15043f2090d5SJed Brown   PetscFunctionBegin;
15050700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
150672ac3e02SJed Brown   if (ts->ops->prestep) {
15073f2090d5SJed Brown     PetscStackPush("TS PreStep function");
15083f2090d5SJed Brown     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
15093f2090d5SJed Brown     PetscStackPop;
1510312ce896SJed Brown   }
15113f2090d5SJed Brown   PetscFunctionReturn(0);
15123f2090d5SJed Brown }
15133f2090d5SJed Brown 
15143f2090d5SJed Brown #undef __FUNCT__
1515e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1516ac226902SBarry Smith /*@C
1517000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
15183f2090d5SJed Brown   called once at the end of each time step.
1519000e7ae3SMatthew Knepley 
15203f9fe445SBarry Smith   Logically Collective on TS
1521000e7ae3SMatthew Knepley 
1522000e7ae3SMatthew Knepley   Input Parameters:
1523000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1524000e7ae3SMatthew Knepley - func - The function
1525000e7ae3SMatthew Knepley 
1526000e7ae3SMatthew Knepley   Calling sequence of func:
1527000e7ae3SMatthew Knepley . func (TS ts);
1528000e7ae3SMatthew Knepley 
1529000e7ae3SMatthew Knepley   Level: intermediate
1530000e7ae3SMatthew Knepley 
1531000e7ae3SMatthew Knepley .keywords: TS, timestep
1532000e7ae3SMatthew Knepley @*/
15337087cfbeSBarry Smith PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1534000e7ae3SMatthew Knepley {
1535000e7ae3SMatthew Knepley   PetscFunctionBegin;
15360700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1537000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1538000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1539000e7ae3SMatthew Knepley }
1540000e7ae3SMatthew Knepley 
1541e74ef692SMatthew Knepley #undef __FUNCT__
15423f2090d5SJed Brown #define __FUNCT__ "TSPostStep"
15433f2090d5SJed Brown /*@C
15443f2090d5SJed Brown   TSPostStep - Runs the user-defined post-step function.
15453f2090d5SJed Brown 
15463f2090d5SJed Brown   Collective on TS
15473f2090d5SJed Brown 
15483f2090d5SJed Brown   Input Parameters:
15493f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
15503f2090d5SJed Brown 
15513f2090d5SJed Brown   Notes:
15523f2090d5SJed Brown   TSPostStep() is typically used within time stepping implementations,
15533f2090d5SJed Brown   so most users would not generally call this routine themselves.
15543f2090d5SJed Brown 
15553f2090d5SJed Brown   Level: developer
15563f2090d5SJed Brown 
15573f2090d5SJed Brown .keywords: TS, timestep
15583f2090d5SJed Brown @*/
15597087cfbeSBarry Smith PetscErrorCode  TSPostStep(TS ts)
15603f2090d5SJed Brown {
15613f2090d5SJed Brown   PetscErrorCode ierr;
15623f2090d5SJed Brown 
15633f2090d5SJed Brown   PetscFunctionBegin;
15640700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
156572ac3e02SJed Brown   if (ts->ops->poststep) {
15663f2090d5SJed Brown     PetscStackPush("TS PostStep function");
15673f2090d5SJed Brown     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
15683f2090d5SJed Brown     PetscStackPop;
156972ac3e02SJed Brown   }
15703f2090d5SJed Brown   PetscFunctionReturn(0);
15713f2090d5SJed Brown }
15723f2090d5SJed Brown 
1573d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1574d763cef2SBarry Smith 
15754a2ae208SSatish Balay #undef __FUNCT__
1576a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet"
1577d763cef2SBarry Smith /*@C
1578a6570f20SBarry Smith    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1579d763cef2SBarry Smith    timestep to display the iteration's  progress.
1580d763cef2SBarry Smith 
15813f9fe445SBarry Smith    Logically Collective on TS
1582d763cef2SBarry Smith 
1583d763cef2SBarry Smith    Input Parameters:
1584d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1585e213d8f1SJed Brown .  monitor - monitoring routine
1586329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1587b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1588b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1589b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1590d763cef2SBarry Smith 
1591e213d8f1SJed Brown    Calling sequence of monitor:
1592e213d8f1SJed Brown $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1593d763cef2SBarry Smith 
1594d763cef2SBarry Smith +    ts - the TS context
1595d763cef2SBarry Smith .    steps - iteration number
15961f06c33eSBarry Smith .    time - current time
1597d763cef2SBarry Smith .    x - current iterate
1598d763cef2SBarry Smith -    mctx - [optional] monitoring context
1599d763cef2SBarry Smith 
1600d763cef2SBarry Smith    Notes:
1601d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1602d763cef2SBarry Smith    already has been loaded.
1603d763cef2SBarry Smith 
1604025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each TS object
1605025f1a04SBarry Smith 
1606d763cef2SBarry Smith    Level: intermediate
1607d763cef2SBarry Smith 
1608d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1609d763cef2SBarry Smith 
1610a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel()
1611d763cef2SBarry Smith @*/
1612c2efdce3SBarry Smith PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1613d763cef2SBarry Smith {
1614d763cef2SBarry Smith   PetscFunctionBegin;
16150700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
161617186662SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1617d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1618329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1619d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1620d763cef2SBarry Smith   PetscFunctionReturn(0);
1621d763cef2SBarry Smith }
1622d763cef2SBarry Smith 
16234a2ae208SSatish Balay #undef __FUNCT__
1624a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel"
1625d763cef2SBarry Smith /*@C
1626a6570f20SBarry Smith    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1627d763cef2SBarry Smith 
16283f9fe445SBarry Smith    Logically Collective on TS
1629d763cef2SBarry Smith 
1630d763cef2SBarry Smith    Input Parameters:
1631d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1632d763cef2SBarry Smith 
1633d763cef2SBarry Smith    Notes:
1634d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1635d763cef2SBarry Smith 
1636d763cef2SBarry Smith    Level: intermediate
1637d763cef2SBarry Smith 
1638d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1639d763cef2SBarry Smith 
1640a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet()
1641d763cef2SBarry Smith @*/
16427087cfbeSBarry Smith PetscErrorCode  TSMonitorCancel(TS ts)
1643d763cef2SBarry Smith {
1644d952e501SBarry Smith   PetscErrorCode ierr;
1645d952e501SBarry Smith   PetscInt       i;
1646d952e501SBarry Smith 
1647d763cef2SBarry Smith   PetscFunctionBegin;
16480700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1649d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1650d952e501SBarry Smith     if (ts->mdestroy[i]) {
16513c4aec1bSBarry Smith       ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
1652d952e501SBarry Smith     }
1653d952e501SBarry Smith   }
1654d763cef2SBarry Smith   ts->numbermonitors = 0;
1655d763cef2SBarry Smith   PetscFunctionReturn(0);
1656d763cef2SBarry Smith }
1657d763cef2SBarry Smith 
16584a2ae208SSatish Balay #undef __FUNCT__
1659a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault"
1660d8e5e3e6SSatish Balay /*@
1661a6570f20SBarry Smith    TSMonitorDefault - Sets the Default monitor
16625516499fSSatish Balay 
16635516499fSSatish Balay    Level: intermediate
166441251cbbSSatish Balay 
16655516499fSSatish Balay .keywords: TS, set, monitor
16665516499fSSatish Balay 
166741251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet()
166841251cbbSSatish Balay @*/
1669649052a6SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1670d763cef2SBarry Smith {
1671dfbe8321SBarry Smith   PetscErrorCode ierr;
1672649052a6SBarry Smith   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1673d132466eSBarry Smith 
1674d763cef2SBarry Smith   PetscFunctionBegin;
1675649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1676971832b6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr);
1677649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1678d763cef2SBarry Smith   PetscFunctionReturn(0);
1679d763cef2SBarry Smith }
1680d763cef2SBarry Smith 
16814a2ae208SSatish Balay #undef __FUNCT__
1682cd652676SJed Brown #define __FUNCT__ "TSSetRetainStages"
1683cd652676SJed Brown /*@
1684cd652676SJed Brown    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
1685cd652676SJed Brown 
1686cd652676SJed Brown    Logically Collective on TS
1687cd652676SJed Brown 
1688cd652676SJed Brown    Input Argument:
1689cd652676SJed Brown .  ts - time stepping context
1690cd652676SJed Brown 
1691cd652676SJed Brown    Output Argument:
1692cd652676SJed Brown .  flg - PETSC_TRUE or PETSC_FALSE
1693cd652676SJed Brown 
1694cd652676SJed Brown    Level: intermediate
1695cd652676SJed Brown 
1696cd652676SJed Brown .keywords: TS, set
1697cd652676SJed Brown 
1698cd652676SJed Brown .seealso: TSInterpolate(), TSSetPostStep()
1699cd652676SJed Brown @*/
1700cd652676SJed Brown PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
1701cd652676SJed Brown {
1702cd652676SJed Brown 
1703cd652676SJed Brown   PetscFunctionBegin;
1704cd652676SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1705cd652676SJed Brown   ts->retain_stages = flg;
1706cd652676SJed Brown   PetscFunctionReturn(0);
1707cd652676SJed Brown }
1708cd652676SJed Brown 
1709cd652676SJed Brown #undef __FUNCT__
1710cd652676SJed Brown #define __FUNCT__ "TSInterpolate"
1711cd652676SJed Brown /*@
1712cd652676SJed Brown    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
1713cd652676SJed Brown 
1714cd652676SJed Brown    Collective on TS
1715cd652676SJed Brown 
1716cd652676SJed Brown    Input Argument:
1717cd652676SJed Brown +  ts - time stepping context
1718cd652676SJed Brown -  t - time to interpolate to
1719cd652676SJed Brown 
1720cd652676SJed Brown    Output Argument:
1721cd652676SJed Brown .  X - state at given time
1722cd652676SJed Brown 
1723cd652676SJed Brown    Notes:
1724cd652676SJed Brown    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
1725cd652676SJed Brown 
1726cd652676SJed Brown    Level: intermediate
1727cd652676SJed Brown 
1728cd652676SJed Brown    Developer Notes:
1729cd652676SJed Brown    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
1730cd652676SJed Brown 
1731cd652676SJed Brown .keywords: TS, set
1732cd652676SJed Brown 
1733cd652676SJed Brown .seealso: TSSetRetainStages(), TSSetPostStep()
1734cd652676SJed Brown @*/
1735cd652676SJed Brown PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
1736cd652676SJed Brown {
1737cd652676SJed Brown   PetscErrorCode ierr;
1738cd652676SJed Brown 
1739cd652676SJed Brown   PetscFunctionBegin;
1740cd652676SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1741d8cd7023SBarry Smith   if (t < ts->ptime - ts->time_step_prev || t > ts->ptime) SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested time %G not in last time steps [%G,%G]",t,ts->ptime-ts->time_step_prev,ts->ptime);
1742cd652676SJed Brown   if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
1743cd652676SJed Brown   ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr);
1744cd652676SJed Brown   PetscFunctionReturn(0);
1745cd652676SJed Brown }
1746cd652676SJed Brown 
1747cd652676SJed Brown #undef __FUNCT__
17484a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1749d763cef2SBarry Smith /*@
17506d9e5789SSean Farley    TSStep - Steps one time step
1751d763cef2SBarry Smith 
1752d763cef2SBarry Smith    Collective on TS
1753d763cef2SBarry Smith 
1754d763cef2SBarry Smith    Input Parameter:
1755d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1756d763cef2SBarry Smith 
17576d9e5789SSean Farley    Level: intermediate
1758d763cef2SBarry Smith 
1759d763cef2SBarry Smith .keywords: TS, timestep, solve
1760d763cef2SBarry Smith 
17616d9e5789SSean Farley .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve()
1762d763cef2SBarry Smith @*/
1763193ac0bcSJed Brown PetscErrorCode  TSStep(TS ts)
1764d763cef2SBarry Smith {
1765362cd11cSLisandro Dalcin   PetscReal      ptime_prev;
1766dfbe8321SBarry Smith   PetscErrorCode ierr;
1767d763cef2SBarry Smith 
1768d763cef2SBarry Smith   PetscFunctionBegin;
17690700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1770d405a339SMatthew Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
1771d405a339SMatthew Knepley 
1772362cd11cSLisandro Dalcin   ts->reason = TS_CONVERGED_ITERATING;
1773362cd11cSLisandro Dalcin 
1774362cd11cSLisandro Dalcin   ptime_prev = ts->ptime;
1775d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1776193ac0bcSJed Brown   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
1777d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1778362cd11cSLisandro Dalcin   ts->time_step_prev = ts->ptime - ptime_prev;
1779362cd11cSLisandro Dalcin 
1780362cd11cSLisandro Dalcin   if (ts->reason < 0) {
1781362cd11cSLisandro Dalcin     if (ts->errorifstepfailed) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed");
1782362cd11cSLisandro Dalcin   } else if (!ts->reason) {
1783362cd11cSLisandro Dalcin     if (ts->steps >= ts->max_steps)
1784362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_ITS;
1785362cd11cSLisandro Dalcin     else if (ts->ptime >= ts->max_time)
1786362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_TIME;
1787362cd11cSLisandro Dalcin   }
1788362cd11cSLisandro Dalcin 
1789d763cef2SBarry Smith   PetscFunctionReturn(0);
1790d763cef2SBarry Smith }
1791d763cef2SBarry Smith 
17924a2ae208SSatish Balay #undef __FUNCT__
17936a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve"
17946a4d4014SLisandro Dalcin /*@
17956a4d4014SLisandro Dalcin    TSSolve - Steps the requested number of timesteps.
17966a4d4014SLisandro Dalcin 
17976a4d4014SLisandro Dalcin    Collective on TS
17986a4d4014SLisandro Dalcin 
17996a4d4014SLisandro Dalcin    Input Parameter:
18006a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
1801362cd11cSLisandro Dalcin -  x - the solution vector
18026a4d4014SLisandro Dalcin 
18035a3a76d0SJed Brown    Output Parameter:
18045a3a76d0SJed Brown .  ftime - time of the state vector x upon completion
18055a3a76d0SJed Brown 
18066a4d4014SLisandro Dalcin    Level: beginner
18076a4d4014SLisandro Dalcin 
18085a3a76d0SJed Brown    Notes:
18095a3a76d0SJed Brown    The final time returned by this function may be different from the time of the internally
18105a3a76d0SJed Brown    held state accessible by TSGetSolution() and TSGetTime() because the method may have
18115a3a76d0SJed Brown    stepped over the final time.
18125a3a76d0SJed Brown 
18136a4d4014SLisandro Dalcin .keywords: TS, timestep, solve
18146a4d4014SLisandro Dalcin 
18156a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep()
18166a4d4014SLisandro Dalcin @*/
18175a3a76d0SJed Brown PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime)
18186a4d4014SLisandro Dalcin {
18194d7d938eSLisandro Dalcin   PetscBool      flg;
18204d7d938eSLisandro Dalcin   char           filename[PETSC_MAX_PATH_LEN];
18214d7d938eSLisandro Dalcin   PetscViewer    viewer;
18226a4d4014SLisandro Dalcin   PetscErrorCode ierr;
1823f22f69f0SBarry Smith 
18246a4d4014SLisandro Dalcin   PetscFunctionBegin;
18250700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1826193ac0bcSJed Brown   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
18275a3a76d0SJed Brown   if (ts->exact_final_time) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
18285a3a76d0SJed Brown     if (!ts->vec_sol || x == ts->vec_sol) {
18295a3a76d0SJed Brown       Vec y;
18305a3a76d0SJed Brown       ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
18315a3a76d0SJed Brown       ierr = TSSetSolution(ts,y);CHKERRQ(ierr);
18325a3a76d0SJed Brown       ierr = VecDestroy(&y);CHKERRQ(ierr); /* grant ownership */
18335a3a76d0SJed Brown     }
1834362cd11cSLisandro Dalcin     ierr = VecCopy(x,ts->vec_sol);CHKERRQ(ierr);
18355a3a76d0SJed Brown   } else {
1836193ac0bcSJed Brown     ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
18375a3a76d0SJed Brown   }
1838b5d403baSSean Farley   ierr = TSSetUp(ts);CHKERRQ(ierr);
18396a4d4014SLisandro Dalcin   /* reset time step and iteration counters */
1840193ac0bcSJed Brown   ts->steps = 0;
1841193ac0bcSJed Brown   ts->linear_its = 0;
1842193ac0bcSJed Brown   ts->nonlinear_its = 0;
1843c610991cSLisandro Dalcin   ts->num_snes_failures = 0;
1844c610991cSLisandro Dalcin   ts->reject = 0;
1845193ac0bcSJed Brown   ts->reason = TS_CONVERGED_ITERATING;
1846193ac0bcSJed Brown 
1847193ac0bcSJed Brown   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
1848193ac0bcSJed Brown     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
18495a3a76d0SJed Brown     ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
1850a5b82c69SSean Farley     if (ftime) *ftime = ts->ptime;
1851193ac0bcSJed Brown   } else {
18526a4d4014SLisandro Dalcin     /* steps the requested number of timesteps. */
1853362cd11cSLisandro Dalcin     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1854362cd11cSLisandro Dalcin     if (ts->steps >= ts->max_steps)
1855362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_ITS;
1856362cd11cSLisandro Dalcin     else if (ts->ptime >= ts->max_time)
1857362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_TIME;
1858e1a7a14fSJed Brown     while (!ts->reason) {
1859193ac0bcSJed Brown       ierr = TSPreStep(ts);CHKERRQ(ierr);
1860193ac0bcSJed Brown       ierr = TSStep(ts);CHKERRQ(ierr);
1861193ac0bcSJed Brown       ierr = TSPostStep(ts);CHKERRQ(ierr);
1862193ac0bcSJed Brown       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1863193ac0bcSJed Brown     }
1864362cd11cSLisandro Dalcin     if (ts->exact_final_time && ts->ptime > ts->max_time) {
1865a43b19c4SJed Brown       ierr = TSInterpolate(ts,ts->max_time,x);CHKERRQ(ierr);
18665a3a76d0SJed Brown       if (ftime) *ftime = ts->max_time;
18670574a7fbSJed Brown     } else {
18680574a7fbSJed Brown       ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
18690574a7fbSJed Brown       if (ftime) *ftime = ts->ptime;
18700574a7fbSJed Brown     }
1871193ac0bcSJed Brown   }
18724d7d938eSLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
18734d7d938eSLisandro Dalcin   if (flg && !PetscPreLoadingOn) {
18744d7d938eSLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr);
18754d7d938eSLisandro Dalcin     ierr = TSView(ts,viewer);CHKERRQ(ierr);
18764d7d938eSLisandro Dalcin     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1877193ac0bcSJed Brown   }
18786a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
18796a4d4014SLisandro Dalcin }
18806a4d4014SLisandro Dalcin 
18816a4d4014SLisandro Dalcin #undef __FUNCT__
18824a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
1883d763cef2SBarry Smith /*
1884d763cef2SBarry Smith      Runs the user provided monitor routines, if they exists.
1885d763cef2SBarry Smith */
1886a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
1887d763cef2SBarry Smith {
18886849ba73SBarry Smith   PetscErrorCode ierr;
1889a7cc72afSBarry Smith   PetscInt       i,n = ts->numbermonitors;
1890d763cef2SBarry Smith 
1891d763cef2SBarry Smith   PetscFunctionBegin;
1892d763cef2SBarry Smith   for (i=0; i<n; i++) {
1893142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
1894d763cef2SBarry Smith   }
1895d763cef2SBarry Smith   PetscFunctionReturn(0);
1896d763cef2SBarry Smith }
1897d763cef2SBarry Smith 
1898d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
1899d763cef2SBarry Smith 
19004a2ae208SSatish Balay #undef __FUNCT__
1901a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate"
1902d763cef2SBarry Smith /*@C
1903a6570f20SBarry Smith    TSMonitorLGCreate - Creates a line graph context for use with
1904d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
1905d763cef2SBarry Smith 
1906d763cef2SBarry Smith    Collective on TS
1907d763cef2SBarry Smith 
1908d763cef2SBarry Smith    Input Parameters:
1909d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
1910d763cef2SBarry Smith .  label - the title to put in the title bar
19117c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
1912d763cef2SBarry Smith -  m, n - the screen width and height in pixels
1913d763cef2SBarry Smith 
1914d763cef2SBarry Smith    Output Parameter:
1915d763cef2SBarry Smith .  draw - the drawing context
1916d763cef2SBarry Smith 
1917d763cef2SBarry Smith    Options Database Key:
1918a6570f20SBarry Smith .  -ts_monitor_draw - automatically sets line graph monitor
1919d763cef2SBarry Smith 
1920d763cef2SBarry Smith    Notes:
1921a6570f20SBarry Smith    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
1922d763cef2SBarry Smith 
1923d763cef2SBarry Smith    Level: intermediate
1924d763cef2SBarry Smith 
19257c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
1926d763cef2SBarry Smith 
1927a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet()
19287c922b88SBarry Smith 
1929d763cef2SBarry Smith @*/
19307087cfbeSBarry Smith PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1931d763cef2SBarry Smith {
1932b0a32e0cSBarry Smith   PetscDraw      win;
1933dfbe8321SBarry Smith   PetscErrorCode ierr;
1934d763cef2SBarry Smith 
1935d763cef2SBarry Smith   PetscFunctionBegin;
1936b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
1937b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
1938b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
1939b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
1940d763cef2SBarry Smith 
194152e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
1942d763cef2SBarry Smith   PetscFunctionReturn(0);
1943d763cef2SBarry Smith }
1944d763cef2SBarry Smith 
19454a2ae208SSatish Balay #undef __FUNCT__
1946a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG"
1947a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
1948d763cef2SBarry Smith {
1949b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
195087828ca2SBarry Smith   PetscReal      x,y = ptime;
1951dfbe8321SBarry Smith   PetscErrorCode ierr;
1952d763cef2SBarry Smith 
1953d763cef2SBarry Smith   PetscFunctionBegin;
19547c922b88SBarry Smith   if (!monctx) {
19557c922b88SBarry Smith     MPI_Comm    comm;
1956b0a32e0cSBarry Smith     PetscViewer viewer;
19577c922b88SBarry Smith 
19587c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
1959b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
1960b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
19617c922b88SBarry Smith   }
19627c922b88SBarry Smith 
1963b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
196487828ca2SBarry Smith   x = (PetscReal)n;
1965b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1966d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
1967b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1968d763cef2SBarry Smith   }
1969d763cef2SBarry Smith   PetscFunctionReturn(0);
1970d763cef2SBarry Smith }
1971d763cef2SBarry Smith 
19724a2ae208SSatish Balay #undef __FUNCT__
1973a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy"
1974d763cef2SBarry Smith /*@C
1975a6570f20SBarry Smith    TSMonitorLGDestroy - Destroys a line graph context that was created
1976a6570f20SBarry Smith    with TSMonitorLGCreate().
1977d763cef2SBarry Smith 
1978b0a32e0cSBarry Smith    Collective on PetscDrawLG
1979d763cef2SBarry Smith 
1980d763cef2SBarry Smith    Input Parameter:
1981d763cef2SBarry Smith .  draw - the drawing context
1982d763cef2SBarry Smith 
1983d763cef2SBarry Smith    Level: intermediate
1984d763cef2SBarry Smith 
1985d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
1986d763cef2SBarry Smith 
1987a6570f20SBarry Smith .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
1988d763cef2SBarry Smith @*/
19896bf464f9SBarry Smith PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
1990d763cef2SBarry Smith {
1991b0a32e0cSBarry Smith   PetscDraw      draw;
1992dfbe8321SBarry Smith   PetscErrorCode ierr;
1993d763cef2SBarry Smith 
1994d763cef2SBarry Smith   PetscFunctionBegin;
19956bf464f9SBarry Smith   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
19966bf464f9SBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
1997b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
1998d763cef2SBarry Smith   PetscFunctionReturn(0);
1999d763cef2SBarry Smith }
2000d763cef2SBarry Smith 
20014a2ae208SSatish Balay #undef __FUNCT__
20024a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
2003d763cef2SBarry Smith /*@
2004d763cef2SBarry Smith    TSGetTime - Gets the current time.
2005d763cef2SBarry Smith 
2006d763cef2SBarry Smith    Not Collective
2007d763cef2SBarry Smith 
2008d763cef2SBarry Smith    Input Parameter:
2009d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
2010d763cef2SBarry Smith 
2011d763cef2SBarry Smith    Output Parameter:
2012d763cef2SBarry Smith .  t  - the current time
2013d763cef2SBarry Smith 
2014d763cef2SBarry Smith    Level: beginner
2015d763cef2SBarry Smith 
2016d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
2017d763cef2SBarry Smith 
2018d763cef2SBarry Smith .keywords: TS, get, time
2019d763cef2SBarry Smith @*/
20207087cfbeSBarry Smith PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
2021d763cef2SBarry Smith {
2022d763cef2SBarry Smith   PetscFunctionBegin;
20230700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
20244482741eSBarry Smith   PetscValidDoublePointer(t,2);
2025d763cef2SBarry Smith   *t = ts->ptime;
2026d763cef2SBarry Smith   PetscFunctionReturn(0);
2027d763cef2SBarry Smith }
2028d763cef2SBarry Smith 
20294a2ae208SSatish Balay #undef __FUNCT__
20306a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime"
20316a4d4014SLisandro Dalcin /*@
20326a4d4014SLisandro Dalcin    TSSetTime - Allows one to reset the time.
20336a4d4014SLisandro Dalcin 
20343f9fe445SBarry Smith    Logically Collective on TS
20356a4d4014SLisandro Dalcin 
20366a4d4014SLisandro Dalcin    Input Parameters:
20376a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
20386a4d4014SLisandro Dalcin -  time - the time
20396a4d4014SLisandro Dalcin 
20406a4d4014SLisandro Dalcin    Level: intermediate
20416a4d4014SLisandro Dalcin 
20426a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration()
20436a4d4014SLisandro Dalcin 
20446a4d4014SLisandro Dalcin .keywords: TS, set, time
20456a4d4014SLisandro Dalcin @*/
20467087cfbeSBarry Smith PetscErrorCode  TSSetTime(TS ts, PetscReal t)
20476a4d4014SLisandro Dalcin {
20486a4d4014SLisandro Dalcin   PetscFunctionBegin;
20490700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2050c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,t,2);
20516a4d4014SLisandro Dalcin   ts->ptime = t;
20526a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
20536a4d4014SLisandro Dalcin }
20546a4d4014SLisandro Dalcin 
20556a4d4014SLisandro Dalcin #undef __FUNCT__
20564a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
2057d763cef2SBarry Smith /*@C
2058d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
2059d763cef2SBarry Smith    TS options in the database.
2060d763cef2SBarry Smith 
20613f9fe445SBarry Smith    Logically Collective on TS
2062d763cef2SBarry Smith 
2063d763cef2SBarry Smith    Input Parameter:
2064d763cef2SBarry Smith +  ts     - The TS context
2065d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
2066d763cef2SBarry Smith 
2067d763cef2SBarry Smith    Notes:
2068d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
2069d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
2070d763cef2SBarry Smith    hyphen.
2071d763cef2SBarry Smith 
2072d763cef2SBarry Smith    Level: advanced
2073d763cef2SBarry Smith 
2074d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
2075d763cef2SBarry Smith 
2076d763cef2SBarry Smith .seealso: TSSetFromOptions()
2077d763cef2SBarry Smith 
2078d763cef2SBarry Smith @*/
20797087cfbeSBarry Smith PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2080d763cef2SBarry Smith {
2081dfbe8321SBarry Smith   PetscErrorCode ierr;
2082089b2837SJed Brown   SNES           snes;
2083d763cef2SBarry Smith 
2084d763cef2SBarry Smith   PetscFunctionBegin;
20850700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2086d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2087089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2088089b2837SJed Brown   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2089d763cef2SBarry Smith   PetscFunctionReturn(0);
2090d763cef2SBarry Smith }
2091d763cef2SBarry Smith 
2092d763cef2SBarry Smith 
20934a2ae208SSatish Balay #undef __FUNCT__
20944a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
2095d763cef2SBarry Smith /*@C
2096d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2097d763cef2SBarry Smith    TS options in the database.
2098d763cef2SBarry Smith 
20993f9fe445SBarry Smith    Logically Collective on TS
2100d763cef2SBarry Smith 
2101d763cef2SBarry Smith    Input Parameter:
2102d763cef2SBarry Smith +  ts     - The TS context
2103d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
2104d763cef2SBarry Smith 
2105d763cef2SBarry Smith    Notes:
2106d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
2107d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
2108d763cef2SBarry Smith    hyphen.
2109d763cef2SBarry Smith 
2110d763cef2SBarry Smith    Level: advanced
2111d763cef2SBarry Smith 
2112d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
2113d763cef2SBarry Smith 
2114d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
2115d763cef2SBarry Smith 
2116d763cef2SBarry Smith @*/
21177087cfbeSBarry Smith PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2118d763cef2SBarry Smith {
2119dfbe8321SBarry Smith   PetscErrorCode ierr;
2120089b2837SJed Brown   SNES           snes;
2121d763cef2SBarry Smith 
2122d763cef2SBarry Smith   PetscFunctionBegin;
21230700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2124d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2125089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2126089b2837SJed Brown   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2127d763cef2SBarry Smith   PetscFunctionReturn(0);
2128d763cef2SBarry Smith }
2129d763cef2SBarry Smith 
21304a2ae208SSatish Balay #undef __FUNCT__
21314a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
2132d763cef2SBarry Smith /*@C
2133d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
2134d763cef2SBarry Smith    TS options in the database.
2135d763cef2SBarry Smith 
2136d763cef2SBarry Smith    Not Collective
2137d763cef2SBarry Smith 
2138d763cef2SBarry Smith    Input Parameter:
2139d763cef2SBarry Smith .  ts - The TS context
2140d763cef2SBarry Smith 
2141d763cef2SBarry Smith    Output Parameter:
2142d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
2143d763cef2SBarry Smith 
2144d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
2145d763cef2SBarry Smith    sufficient length to hold the prefix.
2146d763cef2SBarry Smith 
2147d763cef2SBarry Smith    Level: intermediate
2148d763cef2SBarry Smith 
2149d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
2150d763cef2SBarry Smith 
2151d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
2152d763cef2SBarry Smith @*/
21537087cfbeSBarry Smith PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2154d763cef2SBarry Smith {
2155dfbe8321SBarry Smith   PetscErrorCode ierr;
2156d763cef2SBarry Smith 
2157d763cef2SBarry Smith   PetscFunctionBegin;
21580700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
21594482741eSBarry Smith   PetscValidPointer(prefix,2);
2160d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2161d763cef2SBarry Smith   PetscFunctionReturn(0);
2162d763cef2SBarry Smith }
2163d763cef2SBarry Smith 
21644a2ae208SSatish Balay #undef __FUNCT__
21654a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
2166d763cef2SBarry Smith /*@C
2167d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2168d763cef2SBarry Smith 
2169d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
2170d763cef2SBarry Smith 
2171d763cef2SBarry Smith    Input Parameter:
2172d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
2173d763cef2SBarry Smith 
2174d763cef2SBarry Smith    Output Parameters:
2175d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
2176d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
2177089b2837SJed Brown .  func - Function to compute the Jacobian of the RHS
2178d763cef2SBarry Smith -  ctx - User-defined context for Jacobian evaluation routine
2179d763cef2SBarry Smith 
2180d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2181d763cef2SBarry Smith 
2182d763cef2SBarry Smith    Level: intermediate
2183d763cef2SBarry Smith 
218426d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2185d763cef2SBarry Smith 
2186d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
2187d763cef2SBarry Smith @*/
2188089b2837SJed Brown PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2189d763cef2SBarry Smith {
2190089b2837SJed Brown   PetscErrorCode ierr;
2191089b2837SJed Brown   SNES           snes;
2192089b2837SJed Brown 
2193d763cef2SBarry Smith   PetscFunctionBegin;
2194089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2195089b2837SJed Brown   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
21963b5f76d0SSean Farley   if (func) *func = ts->userops->rhsjacobian;
219726d46c62SHong Zhang   if (ctx) *ctx = ts->jacP;
2198d763cef2SBarry Smith   PetscFunctionReturn(0);
2199d763cef2SBarry Smith }
2200d763cef2SBarry Smith 
22011713a123SBarry Smith #undef __FUNCT__
22022eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian"
22032eca1d9cSJed Brown /*@C
22042eca1d9cSJed Brown    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
22052eca1d9cSJed Brown 
22062eca1d9cSJed Brown    Not Collective, but parallel objects are returned if TS is parallel
22072eca1d9cSJed Brown 
22082eca1d9cSJed Brown    Input Parameter:
22092eca1d9cSJed Brown .  ts  - The TS context obtained from TSCreate()
22102eca1d9cSJed Brown 
22112eca1d9cSJed Brown    Output Parameters:
22122eca1d9cSJed Brown +  A   - The Jacobian of F(t,U,U_t)
22132eca1d9cSJed Brown .  B   - The preconditioner matrix, often the same as A
22142eca1d9cSJed Brown .  f   - The function to compute the matrices
22152eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine
22162eca1d9cSJed Brown 
22172eca1d9cSJed Brown    Notes: You can pass in PETSC_NULL for any return argument you do not need.
22182eca1d9cSJed Brown 
22192eca1d9cSJed Brown    Level: advanced
22202eca1d9cSJed Brown 
22212eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
22222eca1d9cSJed Brown 
22232eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian
22242eca1d9cSJed Brown @*/
22257087cfbeSBarry Smith PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
22262eca1d9cSJed Brown {
2227089b2837SJed Brown   PetscErrorCode ierr;
2228089b2837SJed Brown   SNES           snes;
2229089b2837SJed Brown 
22302eca1d9cSJed Brown   PetscFunctionBegin;
2231089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2232089b2837SJed Brown   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
22333b5f76d0SSean Farley   if (f) *f = ts->userops->ijacobian;
22342eca1d9cSJed Brown   if (ctx) *ctx = ts->jacP;
22352eca1d9cSJed Brown   PetscFunctionReturn(0);
22362eca1d9cSJed Brown }
22372eca1d9cSJed Brown 
22386083293cSBarry Smith typedef struct {
22396083293cSBarry Smith   PetscViewer viewer;
22406083293cSBarry Smith   Vec         initialsolution;
22416083293cSBarry Smith   PetscBool   showinitial;
22426083293cSBarry Smith } TSMonitorSolutionCtx;
22436083293cSBarry Smith 
22442eca1d9cSJed Brown #undef __FUNCT__
2245a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution"
22461713a123SBarry Smith /*@C
2247a6570f20SBarry Smith    TSMonitorSolution - Monitors progress of the TS solvers by calling
22481713a123SBarry Smith    VecView() for the solution at each timestep
22491713a123SBarry Smith 
22501713a123SBarry Smith    Collective on TS
22511713a123SBarry Smith 
22521713a123SBarry Smith    Input Parameters:
22531713a123SBarry Smith +  ts - the TS context
22541713a123SBarry Smith .  step - current time-step
2255142b95e3SSatish Balay .  ptime - current time
22561713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
22571713a123SBarry Smith 
22581713a123SBarry Smith    Level: intermediate
22591713a123SBarry Smith 
22601713a123SBarry Smith .keywords: TS,  vector, monitor, view
22611713a123SBarry Smith 
2262a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
22631713a123SBarry Smith @*/
22647087cfbeSBarry Smith PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
22651713a123SBarry Smith {
2266dfbe8321SBarry Smith   PetscErrorCode       ierr;
22676083293cSBarry Smith   TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy;
22681713a123SBarry Smith 
22691713a123SBarry Smith   PetscFunctionBegin;
22706083293cSBarry Smith   if (!step && ictx->showinitial) {
22716083293cSBarry Smith     if (!ictx->initialsolution) {
22726083293cSBarry Smith       ierr = VecDuplicate(x,&ictx->initialsolution);CHKERRQ(ierr);
22731713a123SBarry Smith     }
22746083293cSBarry Smith     ierr = VecCopy(x,ictx->initialsolution);CHKERRQ(ierr);
22756083293cSBarry Smith   }
22766083293cSBarry Smith   if (ictx->showinitial) {
22776083293cSBarry Smith     PetscReal pause;
22786083293cSBarry Smith     ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr);
22796083293cSBarry Smith     ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr);
22806083293cSBarry Smith     ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr);
22816083293cSBarry Smith     ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr);
22826083293cSBarry Smith     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr);
22836083293cSBarry Smith   }
22846083293cSBarry Smith   ierr = VecView(x,ictx->viewer);CHKERRQ(ierr);
22856083293cSBarry Smith   if (ictx->showinitial) {
22866083293cSBarry Smith     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr);
22876083293cSBarry Smith   }
22881713a123SBarry Smith   PetscFunctionReturn(0);
22891713a123SBarry Smith }
22901713a123SBarry Smith 
22911713a123SBarry Smith 
22926c699258SBarry Smith #undef __FUNCT__
22936083293cSBarry Smith #define __FUNCT__ "TSMonitorSolutionDestroy"
22946083293cSBarry Smith /*@C
22956083293cSBarry Smith    TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution
22966083293cSBarry Smith 
22976083293cSBarry Smith    Collective on TS
22986083293cSBarry Smith 
22996083293cSBarry Smith    Input Parameters:
23006083293cSBarry Smith .    ctx - the monitor context
23016083293cSBarry Smith 
23026083293cSBarry Smith    Level: intermediate
23036083293cSBarry Smith 
23046083293cSBarry Smith .keywords: TS,  vector, monitor, view
23056083293cSBarry Smith 
23066083293cSBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
23076083293cSBarry Smith @*/
23086083293cSBarry Smith PetscErrorCode  TSMonitorSolutionDestroy(void **ctx)
23096083293cSBarry Smith {
23106083293cSBarry Smith   PetscErrorCode       ierr;
23116083293cSBarry Smith   TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx;
23126083293cSBarry Smith 
23136083293cSBarry Smith   PetscFunctionBegin;
23146083293cSBarry Smith   ierr = PetscViewerDestroy(&ictx->viewer);CHKERRQ(ierr);
23156083293cSBarry Smith   ierr = VecDestroy(&ictx->initialsolution);CHKERRQ(ierr);
23166083293cSBarry Smith   ierr = PetscFree(ictx);CHKERRQ(ierr);
23176083293cSBarry Smith   PetscFunctionReturn(0);
23186083293cSBarry Smith }
23196083293cSBarry Smith 
23206083293cSBarry Smith #undef __FUNCT__
23216083293cSBarry Smith #define __FUNCT__ "TSMonitorSolutionCreate"
23226083293cSBarry Smith /*@C
23236083293cSBarry Smith    TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution
23246083293cSBarry Smith 
23256083293cSBarry Smith    Collective on TS
23266083293cSBarry Smith 
23276083293cSBarry Smith    Input Parameter:
23286083293cSBarry Smith .    ts - time-step context
23296083293cSBarry Smith 
23306083293cSBarry Smith    Output Patameter:
23316083293cSBarry Smith .    ctx - the monitor context
23326083293cSBarry Smith 
23336083293cSBarry Smith    Level: intermediate
23346083293cSBarry Smith 
23356083293cSBarry Smith .keywords: TS,  vector, monitor, view
23366083293cSBarry Smith 
23376083293cSBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
23386083293cSBarry Smith @*/
23396083293cSBarry Smith PetscErrorCode  TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx)
23406083293cSBarry Smith {
23416083293cSBarry Smith   PetscErrorCode       ierr;
23426083293cSBarry Smith   TSMonitorSolutionCtx *ictx;
23436083293cSBarry Smith 
23446083293cSBarry Smith   PetscFunctionBegin;
23456083293cSBarry Smith   ierr = PetscNew(TSMonitorSolutionCtx,&ictx);CHKERRQ(ierr);
23466083293cSBarry Smith   *ctx = (void*)ictx;
23476083293cSBarry Smith   if (!viewer) {
23486083293cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
23496083293cSBarry Smith   }
23506083293cSBarry Smith   ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
23516083293cSBarry Smith   ictx->viewer      = viewer;
23526083293cSBarry Smith   ictx->showinitial = PETSC_FALSE;
23536083293cSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);CHKERRQ(ierr);
23546083293cSBarry Smith   PetscFunctionReturn(0);
23556083293cSBarry Smith }
23566083293cSBarry Smith 
23576083293cSBarry Smith #undef __FUNCT__
23586c699258SBarry Smith #define __FUNCT__ "TSSetDM"
23596c699258SBarry Smith /*@
23606c699258SBarry Smith    TSSetDM - Sets the DM that may be used by some preconditioners
23616c699258SBarry Smith 
23623f9fe445SBarry Smith    Logically Collective on TS and DM
23636c699258SBarry Smith 
23646c699258SBarry Smith    Input Parameters:
23656c699258SBarry Smith +  ts - the preconditioner context
23666c699258SBarry Smith -  dm - the dm
23676c699258SBarry Smith 
23686c699258SBarry Smith    Level: intermediate
23696c699258SBarry Smith 
23706c699258SBarry Smith 
23716c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
23726c699258SBarry Smith @*/
23737087cfbeSBarry Smith PetscErrorCode  TSSetDM(TS ts,DM dm)
23746c699258SBarry Smith {
23756c699258SBarry Smith   PetscErrorCode ierr;
2376089b2837SJed Brown   SNES           snes;
23776c699258SBarry Smith 
23786c699258SBarry Smith   PetscFunctionBegin;
23790700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
238070663e4aSLisandro Dalcin   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
23816bf464f9SBarry Smith   ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
23826c699258SBarry Smith   ts->dm = dm;
2383089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2384089b2837SJed Brown   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
23856c699258SBarry Smith   PetscFunctionReturn(0);
23866c699258SBarry Smith }
23876c699258SBarry Smith 
23886c699258SBarry Smith #undef __FUNCT__
23896c699258SBarry Smith #define __FUNCT__ "TSGetDM"
23906c699258SBarry Smith /*@
23916c699258SBarry Smith    TSGetDM - Gets the DM that may be used by some preconditioners
23926c699258SBarry Smith 
23933f9fe445SBarry Smith    Not Collective
23946c699258SBarry Smith 
23956c699258SBarry Smith    Input Parameter:
23966c699258SBarry Smith . ts - the preconditioner context
23976c699258SBarry Smith 
23986c699258SBarry Smith    Output Parameter:
23996c699258SBarry Smith .  dm - the dm
24006c699258SBarry Smith 
24016c699258SBarry Smith    Level: intermediate
24026c699258SBarry Smith 
24036c699258SBarry Smith 
24046c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
24056c699258SBarry Smith @*/
24067087cfbeSBarry Smith PetscErrorCode  TSGetDM(TS ts,DM *dm)
24076c699258SBarry Smith {
24086c699258SBarry Smith   PetscFunctionBegin;
24090700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
24106c699258SBarry Smith   *dm = ts->dm;
24116c699258SBarry Smith   PetscFunctionReturn(0);
24126c699258SBarry Smith }
24131713a123SBarry Smith 
24140f5c6efeSJed Brown #undef __FUNCT__
24150f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction"
24160f5c6efeSJed Brown /*@
24170f5c6efeSJed Brown    SNESTSFormFunction - Function to evaluate nonlinear residual
24180f5c6efeSJed Brown 
24193f9fe445SBarry Smith    Logically Collective on SNES
24200f5c6efeSJed Brown 
24210f5c6efeSJed Brown    Input Parameter:
2422d42a1c89SJed Brown + snes - nonlinear solver
24230f5c6efeSJed Brown . X - the current state at which to evaluate the residual
2424d42a1c89SJed Brown - ctx - user context, must be a TS
24250f5c6efeSJed Brown 
24260f5c6efeSJed Brown    Output Parameter:
24270f5c6efeSJed Brown . F - the nonlinear residual
24280f5c6efeSJed Brown 
24290f5c6efeSJed Brown    Notes:
24300f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
24310f5c6efeSJed Brown    It is most frequently passed to MatFDColoringSetFunction().
24320f5c6efeSJed Brown 
24330f5c6efeSJed Brown    Level: advanced
24340f5c6efeSJed Brown 
24350f5c6efeSJed Brown .seealso: SNESSetFunction(), MatFDColoringSetFunction()
24360f5c6efeSJed Brown @*/
24377087cfbeSBarry Smith PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
24380f5c6efeSJed Brown {
24390f5c6efeSJed Brown   TS ts = (TS)ctx;
24400f5c6efeSJed Brown   PetscErrorCode ierr;
24410f5c6efeSJed Brown 
24420f5c6efeSJed Brown   PetscFunctionBegin;
24430f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
24440f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
24450f5c6efeSJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
24460f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
24470f5c6efeSJed Brown   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
24480f5c6efeSJed Brown   PetscFunctionReturn(0);
24490f5c6efeSJed Brown }
24500f5c6efeSJed Brown 
24510f5c6efeSJed Brown #undef __FUNCT__
24520f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian"
24530f5c6efeSJed Brown /*@
24540f5c6efeSJed Brown    SNESTSFormJacobian - Function to evaluate the Jacobian
24550f5c6efeSJed Brown 
24560f5c6efeSJed Brown    Collective on SNES
24570f5c6efeSJed Brown 
24580f5c6efeSJed Brown    Input Parameter:
24590f5c6efeSJed Brown + snes - nonlinear solver
24600f5c6efeSJed Brown . X - the current state at which to evaluate the residual
24610f5c6efeSJed Brown - ctx - user context, must be a TS
24620f5c6efeSJed Brown 
24630f5c6efeSJed Brown    Output Parameter:
24640f5c6efeSJed Brown + A - the Jacobian
24650f5c6efeSJed Brown . B - the preconditioning matrix (may be the same as A)
24660f5c6efeSJed Brown - flag - indicates any structure change in the matrix
24670f5c6efeSJed Brown 
24680f5c6efeSJed Brown    Notes:
24690f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
24700f5c6efeSJed Brown 
24710f5c6efeSJed Brown    Level: developer
24720f5c6efeSJed Brown 
24730f5c6efeSJed Brown .seealso: SNESSetJacobian()
24740f5c6efeSJed Brown @*/
24757087cfbeSBarry Smith PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
24760f5c6efeSJed Brown {
24770f5c6efeSJed Brown   TS ts = (TS)ctx;
24780f5c6efeSJed Brown   PetscErrorCode ierr;
24790f5c6efeSJed Brown 
24800f5c6efeSJed Brown   PetscFunctionBegin;
24810f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
24820f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
24830f5c6efeSJed Brown   PetscValidPointer(A,3);
24840f5c6efeSJed Brown   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
24850f5c6efeSJed Brown   PetscValidPointer(B,4);
24860f5c6efeSJed Brown   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
24870f5c6efeSJed Brown   PetscValidPointer(flag,5);
24880f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
24890f5c6efeSJed Brown   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
24900f5c6efeSJed Brown   PetscFunctionReturn(0);
24910f5c6efeSJed Brown }
2492325fc9f4SBarry Smith 
24930e4ef248SJed Brown #undef __FUNCT__
24940e4ef248SJed Brown #define __FUNCT__ "TSComputeRHSFunctionLinear"
24950e4ef248SJed Brown /*@C
24960e4ef248SJed Brown    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only
24970e4ef248SJed Brown 
24980e4ef248SJed Brown    Collective on TS
24990e4ef248SJed Brown 
25000e4ef248SJed Brown    Input Arguments:
25010e4ef248SJed Brown +  ts - time stepping context
25020e4ef248SJed Brown .  t - time at which to evaluate
25030e4ef248SJed Brown .  X - state at which to evaluate
25040e4ef248SJed Brown -  ctx - context
25050e4ef248SJed Brown 
25060e4ef248SJed Brown    Output Arguments:
25070e4ef248SJed Brown .  F - right hand side
25080e4ef248SJed Brown 
25090e4ef248SJed Brown    Level: intermediate
25100e4ef248SJed Brown 
25110e4ef248SJed Brown    Notes:
25120e4ef248SJed Brown    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
25130e4ef248SJed Brown    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
25140e4ef248SJed Brown 
25150e4ef248SJed Brown .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
25160e4ef248SJed Brown @*/
25170e4ef248SJed Brown PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
25180e4ef248SJed Brown {
25190e4ef248SJed Brown   PetscErrorCode ierr;
25200e4ef248SJed Brown   Mat Arhs,Brhs;
25210e4ef248SJed Brown   MatStructure flg2;
25220e4ef248SJed Brown 
25230e4ef248SJed Brown   PetscFunctionBegin;
25240e4ef248SJed Brown   ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
25250e4ef248SJed Brown   ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
25260e4ef248SJed Brown   ierr = MatMult(Arhs,X,F);CHKERRQ(ierr);
25270e4ef248SJed Brown   PetscFunctionReturn(0);
25280e4ef248SJed Brown }
25290e4ef248SJed Brown 
25300e4ef248SJed Brown #undef __FUNCT__
25310e4ef248SJed Brown #define __FUNCT__ "TSComputeRHSJacobianConstant"
25320e4ef248SJed Brown /*@C
25330e4ef248SJed Brown    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
25340e4ef248SJed Brown 
25350e4ef248SJed Brown    Collective on TS
25360e4ef248SJed Brown 
25370e4ef248SJed Brown    Input Arguments:
25380e4ef248SJed Brown +  ts - time stepping context
25390e4ef248SJed Brown .  t - time at which to evaluate
25400e4ef248SJed Brown .  X - state at which to evaluate
25410e4ef248SJed Brown -  ctx - context
25420e4ef248SJed Brown 
25430e4ef248SJed Brown    Output Arguments:
25440e4ef248SJed Brown +  A - pointer to operator
25450e4ef248SJed Brown .  B - pointer to preconditioning matrix
25460e4ef248SJed Brown -  flg - matrix structure flag
25470e4ef248SJed Brown 
25480e4ef248SJed Brown    Level: intermediate
25490e4ef248SJed Brown 
25500e4ef248SJed Brown    Notes:
25510e4ef248SJed Brown    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
25520e4ef248SJed Brown 
25530e4ef248SJed Brown .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
25540e4ef248SJed Brown @*/
25550e4ef248SJed Brown PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
25560e4ef248SJed Brown {
25570e4ef248SJed Brown 
25580e4ef248SJed Brown   PetscFunctionBegin;
25590e4ef248SJed Brown   *flg = SAME_PRECONDITIONER;
25600e4ef248SJed Brown   PetscFunctionReturn(0);
25610e4ef248SJed Brown }
25620e4ef248SJed Brown 
25630026cea9SSean Farley #undef __FUNCT__
25640026cea9SSean Farley #define __FUNCT__ "TSComputeIFunctionLinear"
25650026cea9SSean Farley /*@C
25660026cea9SSean Farley    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
25670026cea9SSean Farley 
25680026cea9SSean Farley    Collective on TS
25690026cea9SSean Farley 
25700026cea9SSean Farley    Input Arguments:
25710026cea9SSean Farley +  ts - time stepping context
25720026cea9SSean Farley .  t - time at which to evaluate
25730026cea9SSean Farley .  X - state at which to evaluate
25740026cea9SSean Farley .  Xdot - time derivative of state vector
25750026cea9SSean Farley -  ctx - context
25760026cea9SSean Farley 
25770026cea9SSean Farley    Output Arguments:
25780026cea9SSean Farley .  F - left hand side
25790026cea9SSean Farley 
25800026cea9SSean Farley    Level: intermediate
25810026cea9SSean Farley 
25820026cea9SSean Farley    Notes:
25830026cea9SSean Farley    The assumption here is that the left hand side is of the form A*Xdot (and not A*Xdot + B*X). For other cases, the
25840026cea9SSean Farley    user is required to write their own TSComputeIFunction.
25850026cea9SSean Farley    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
25860026cea9SSean Farley    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
25870026cea9SSean Farley 
25880026cea9SSean Farley .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
25890026cea9SSean Farley @*/
25900026cea9SSean Farley PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
25910026cea9SSean Farley {
25920026cea9SSean Farley   PetscErrorCode ierr;
25930026cea9SSean Farley   Mat A,B;
25940026cea9SSean Farley   MatStructure flg2;
25950026cea9SSean Farley 
25960026cea9SSean Farley   PetscFunctionBegin;
25970026cea9SSean Farley   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
25980026cea9SSean Farley   ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr);
25990026cea9SSean Farley   ierr = MatMult(A,Xdot,F);CHKERRQ(ierr);
26000026cea9SSean Farley   PetscFunctionReturn(0);
26010026cea9SSean Farley }
26020026cea9SSean Farley 
26030026cea9SSean Farley #undef __FUNCT__
26040026cea9SSean Farley #define __FUNCT__ "TSComputeIJacobianConstant"
26050026cea9SSean Farley /*@C
26060026cea9SSean Farley    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
26070026cea9SSean Farley 
26080026cea9SSean Farley    Collective on TS
26090026cea9SSean Farley 
26100026cea9SSean Farley    Input Arguments:
26110026cea9SSean Farley +  ts - time stepping context
26120026cea9SSean Farley .  t - time at which to evaluate
26130026cea9SSean Farley .  X - state at which to evaluate
26140026cea9SSean Farley .  Xdot - time derivative of state vector
26150026cea9SSean Farley .  shift - shift to apply
26160026cea9SSean Farley -  ctx - context
26170026cea9SSean Farley 
26180026cea9SSean Farley    Output Arguments:
26190026cea9SSean Farley +  A - pointer to operator
26200026cea9SSean Farley .  B - pointer to preconditioning matrix
26210026cea9SSean Farley -  flg - matrix structure flag
26220026cea9SSean Farley 
26230026cea9SSean Farley    Level: intermediate
26240026cea9SSean Farley 
26250026cea9SSean Farley    Notes:
26260026cea9SSean Farley    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
26270026cea9SSean Farley 
26280026cea9SSean Farley .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
26290026cea9SSean Farley @*/
26300026cea9SSean Farley PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
26310026cea9SSean Farley {
26320026cea9SSean Farley 
26330026cea9SSean Farley   PetscFunctionBegin;
26340026cea9SSean Farley   *flg = SAME_PRECONDITIONER;
26350026cea9SSean Farley   PetscFunctionReturn(0);
26360026cea9SSean Farley }
26370026cea9SSean Farley 
26384af1b03aSJed Brown 
26394af1b03aSJed Brown #undef __FUNCT__
26404af1b03aSJed Brown #define __FUNCT__ "TSGetConvergedReason"
26414af1b03aSJed Brown /*@
26424af1b03aSJed Brown    TSGetConvergedReason - Gets the reason the TS iteration was stopped.
26434af1b03aSJed Brown 
26444af1b03aSJed Brown    Not Collective
26454af1b03aSJed Brown 
26464af1b03aSJed Brown    Input Parameter:
26474af1b03aSJed Brown .  ts - the TS context
26484af1b03aSJed Brown 
26494af1b03aSJed Brown    Output Parameter:
26504af1b03aSJed Brown .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
26514af1b03aSJed Brown             manual pages for the individual convergence tests for complete lists
26524af1b03aSJed Brown 
26534af1b03aSJed Brown    Level: intermediate
26544af1b03aSJed Brown 
2655cd652676SJed Brown    Notes:
2656cd652676SJed Brown    Can only be called after the call to TSSolve() is complete.
26574af1b03aSJed Brown 
26584af1b03aSJed Brown .keywords: TS, nonlinear, set, convergence, test
26594af1b03aSJed Brown 
26604af1b03aSJed Brown .seealso: TSSetConvergenceTest(), TSConvergedReason
26614af1b03aSJed Brown @*/
26624af1b03aSJed Brown PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
26634af1b03aSJed Brown {
26644af1b03aSJed Brown   PetscFunctionBegin;
26654af1b03aSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
26664af1b03aSJed Brown   PetscValidPointer(reason,2);
26674af1b03aSJed Brown   *reason = ts->reason;
26684af1b03aSJed Brown   PetscFunctionReturn(0);
26694af1b03aSJed Brown }
26704af1b03aSJed Brown 
2671d6ebe24aSShri Abhyankar 
2672d6ebe24aSShri Abhyankar #undef __FUNCT__
2673d6ebe24aSShri Abhyankar #define __FUNCT__ "TSVISetVariableBounds"
2674d6ebe24aSShri Abhyankar /*@
2675d6ebe24aSShri Abhyankar    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
2676d6ebe24aSShri Abhyankar 
2677d6ebe24aSShri Abhyankar    Input Parameters:
2678d6ebe24aSShri Abhyankar .  ts   - the TS context.
2679d6ebe24aSShri Abhyankar .  xl   - lower bound.
2680d6ebe24aSShri Abhyankar .  xu   - upper bound.
2681d6ebe24aSShri Abhyankar 
2682d6ebe24aSShri Abhyankar    Notes:
2683d6ebe24aSShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
2684d6ebe24aSShri Abhyankar    SNES_VI_INF and SNES_VI_NINF respectively during SNESSetUp().
2685d6ebe24aSShri Abhyankar 
2686*2bd2b0e6SSatish Balay    Level: advanced
2687*2bd2b0e6SSatish Balay 
2688d6ebe24aSShri Abhyankar @*/
2689d6ebe24aSShri Abhyankar PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
2690d6ebe24aSShri Abhyankar {
2691d6ebe24aSShri Abhyankar   PetscErrorCode ierr;
2692d6ebe24aSShri Abhyankar   SNES           snes;
2693d6ebe24aSShri Abhyankar 
2694d6ebe24aSShri Abhyankar   PetscFunctionBegin;
2695d6ebe24aSShri Abhyankar   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2696d6ebe24aSShri Abhyankar   ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr);
2697d6ebe24aSShri Abhyankar   PetscFunctionReturn(0);
2698d6ebe24aSShri Abhyankar }
2699d6ebe24aSShri Abhyankar 
2700325fc9f4SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
2701c6db04a5SJed Brown #include <mex.h>
2702325fc9f4SBarry Smith 
2703325fc9f4SBarry Smith typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
2704325fc9f4SBarry Smith 
2705325fc9f4SBarry Smith #undef __FUNCT__
2706325fc9f4SBarry Smith #define __FUNCT__ "TSComputeFunction_Matlab"
2707325fc9f4SBarry Smith /*
2708325fc9f4SBarry Smith    TSComputeFunction_Matlab - Calls the function that has been set with
2709325fc9f4SBarry Smith                          TSSetFunctionMatlab().
2710325fc9f4SBarry Smith 
2711325fc9f4SBarry Smith    Collective on TS
2712325fc9f4SBarry Smith 
2713325fc9f4SBarry Smith    Input Parameters:
2714325fc9f4SBarry Smith +  snes - the TS context
2715325fc9f4SBarry Smith -  x - input vector
2716325fc9f4SBarry Smith 
2717325fc9f4SBarry Smith    Output Parameter:
2718325fc9f4SBarry Smith .  y - function vector, as set by TSSetFunction()
2719325fc9f4SBarry Smith 
2720325fc9f4SBarry Smith    Notes:
2721325fc9f4SBarry Smith    TSComputeFunction() is typically used within nonlinear solvers
2722325fc9f4SBarry Smith    implementations, so most users would not generally call this routine
2723325fc9f4SBarry Smith    themselves.
2724325fc9f4SBarry Smith 
2725325fc9f4SBarry Smith    Level: developer
2726325fc9f4SBarry Smith 
2727325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
2728325fc9f4SBarry Smith 
2729325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
2730325fc9f4SBarry Smith */
27317087cfbeSBarry Smith PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
2732325fc9f4SBarry Smith {
2733325fc9f4SBarry Smith   PetscErrorCode   ierr;
2734325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2735325fc9f4SBarry Smith   int              nlhs = 1,nrhs = 7;
2736325fc9f4SBarry Smith   mxArray          *plhs[1],*prhs[7];
2737325fc9f4SBarry Smith   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
2738325fc9f4SBarry Smith 
2739325fc9f4SBarry Smith   PetscFunctionBegin;
2740325fc9f4SBarry Smith   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
2741325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2742325fc9f4SBarry Smith   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
2743325fc9f4SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
2744325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,x,3);
2745325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,y,5);
2746325fc9f4SBarry Smith 
2747325fc9f4SBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
2748325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2749880f3077SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
2750325fc9f4SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
2751325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
2752325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar(time);
2753325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
2754325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2755325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)ly);
2756325fc9f4SBarry Smith   prhs[5] =  mxCreateString(sctx->funcname);
2757325fc9f4SBarry Smith   prhs[6] =  sctx->ctx;
2758325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
2759325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2760325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
2761325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
2762325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
2763325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
2764325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
2765325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
2766325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
2767325fc9f4SBarry Smith   PetscFunctionReturn(0);
2768325fc9f4SBarry Smith }
2769325fc9f4SBarry Smith 
2770325fc9f4SBarry Smith 
2771325fc9f4SBarry Smith #undef __FUNCT__
2772325fc9f4SBarry Smith #define __FUNCT__ "TSSetFunctionMatlab"
2773325fc9f4SBarry Smith /*
2774325fc9f4SBarry Smith    TSSetFunctionMatlab - Sets the function evaluation routine and function
2775325fc9f4SBarry Smith    vector for use by the TS routines in solving ODEs
2776e3c5b3baSBarry Smith    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
2777325fc9f4SBarry Smith 
2778325fc9f4SBarry Smith    Logically Collective on TS
2779325fc9f4SBarry Smith 
2780325fc9f4SBarry Smith    Input Parameters:
2781325fc9f4SBarry Smith +  ts - the TS context
2782325fc9f4SBarry Smith -  func - function evaluation routine
2783325fc9f4SBarry Smith 
2784325fc9f4SBarry Smith    Calling sequence of func:
2785325fc9f4SBarry Smith $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
2786325fc9f4SBarry Smith 
2787325fc9f4SBarry Smith    Level: beginner
2788325fc9f4SBarry Smith 
2789325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
2790325fc9f4SBarry Smith 
2791325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2792325fc9f4SBarry Smith */
2793cdcf91faSSean Farley PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
2794325fc9f4SBarry Smith {
2795325fc9f4SBarry Smith   PetscErrorCode  ierr;
2796325fc9f4SBarry Smith   TSMatlabContext *sctx;
2797325fc9f4SBarry Smith 
2798325fc9f4SBarry Smith   PetscFunctionBegin;
2799325fc9f4SBarry Smith   /* currently sctx is memory bleed */
2800325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2801325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2802325fc9f4SBarry Smith   /*
2803325fc9f4SBarry Smith      This should work, but it doesn't
2804325fc9f4SBarry Smith   sctx->ctx = ctx;
2805325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
2806325fc9f4SBarry Smith   */
2807325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
2808cdcf91faSSean Farley   ierr = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
2809325fc9f4SBarry Smith   PetscFunctionReturn(0);
2810325fc9f4SBarry Smith }
2811325fc9f4SBarry Smith 
2812325fc9f4SBarry Smith #undef __FUNCT__
2813325fc9f4SBarry Smith #define __FUNCT__ "TSComputeJacobian_Matlab"
2814325fc9f4SBarry Smith /*
2815325fc9f4SBarry Smith    TSComputeJacobian_Matlab - Calls the function that has been set with
2816325fc9f4SBarry Smith                          TSSetJacobianMatlab().
2817325fc9f4SBarry Smith 
2818325fc9f4SBarry Smith    Collective on TS
2819325fc9f4SBarry Smith 
2820325fc9f4SBarry Smith    Input Parameters:
2821cdcf91faSSean Farley +  ts - the TS context
2822325fc9f4SBarry Smith .  x - input vector
2823325fc9f4SBarry Smith .  A, B - the matrices
2824325fc9f4SBarry Smith -  ctx - user context
2825325fc9f4SBarry Smith 
2826325fc9f4SBarry Smith    Output Parameter:
2827325fc9f4SBarry Smith .  flag - structure of the matrix
2828325fc9f4SBarry Smith 
2829325fc9f4SBarry Smith    Level: developer
2830325fc9f4SBarry Smith 
2831325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
2832325fc9f4SBarry Smith 
2833325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
2834325fc9f4SBarry Smith @*/
2835cdcf91faSSean Farley PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
2836325fc9f4SBarry Smith {
2837325fc9f4SBarry Smith   PetscErrorCode  ierr;
2838325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2839325fc9f4SBarry Smith   int             nlhs = 2,nrhs = 9;
2840325fc9f4SBarry Smith   mxArray         *plhs[2],*prhs[9];
2841325fc9f4SBarry Smith   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
2842325fc9f4SBarry Smith 
2843325fc9f4SBarry Smith   PetscFunctionBegin;
2844cdcf91faSSean Farley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2845325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2846325fc9f4SBarry Smith 
2847325fc9f4SBarry Smith   /* call Matlab function in ctx with arguments x and y */
2848325fc9f4SBarry Smith 
2849cdcf91faSSean Farley   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
2850325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2851325fc9f4SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
2852325fc9f4SBarry Smith   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
2853325fc9f4SBarry Smith   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
2854325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
2855325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)time);
2856325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
2857325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
2858325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)shift);
2859325fc9f4SBarry Smith   prhs[5] =  mxCreateDoubleScalar((double)lA);
2860325fc9f4SBarry Smith   prhs[6] =  mxCreateDoubleScalar((double)lB);
2861325fc9f4SBarry Smith   prhs[7] =  mxCreateString(sctx->funcname);
2862325fc9f4SBarry Smith   prhs[8] =  sctx->ctx;
2863325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
2864325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2865325fc9f4SBarry Smith   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
2866325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
2867325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
2868325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
2869325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
2870325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
2871325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
2872325fc9f4SBarry Smith   mxDestroyArray(prhs[6]);
2873325fc9f4SBarry Smith   mxDestroyArray(prhs[7]);
2874325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
2875325fc9f4SBarry Smith   mxDestroyArray(plhs[1]);
2876325fc9f4SBarry Smith   PetscFunctionReturn(0);
2877325fc9f4SBarry Smith }
2878325fc9f4SBarry Smith 
2879325fc9f4SBarry Smith 
2880325fc9f4SBarry Smith #undef __FUNCT__
2881325fc9f4SBarry Smith #define __FUNCT__ "TSSetJacobianMatlab"
2882325fc9f4SBarry Smith /*
2883325fc9f4SBarry Smith    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
2884e3c5b3baSBarry 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
2885325fc9f4SBarry Smith 
2886325fc9f4SBarry Smith    Logically Collective on TS
2887325fc9f4SBarry Smith 
2888325fc9f4SBarry Smith    Input Parameters:
2889cdcf91faSSean Farley +  ts - the TS context
2890325fc9f4SBarry Smith .  A,B - Jacobian matrices
2891325fc9f4SBarry Smith .  func - function evaluation routine
2892325fc9f4SBarry Smith -  ctx - user context
2893325fc9f4SBarry Smith 
2894325fc9f4SBarry Smith    Calling sequence of func:
2895cdcf91faSSean Farley $    flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
2896325fc9f4SBarry Smith 
2897325fc9f4SBarry Smith 
2898325fc9f4SBarry Smith    Level: developer
2899325fc9f4SBarry Smith 
2900325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
2901325fc9f4SBarry Smith 
2902325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2903325fc9f4SBarry Smith */
2904cdcf91faSSean Farley PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
2905325fc9f4SBarry Smith {
2906325fc9f4SBarry Smith   PetscErrorCode    ierr;
2907325fc9f4SBarry Smith   TSMatlabContext *sctx;
2908325fc9f4SBarry Smith 
2909325fc9f4SBarry Smith   PetscFunctionBegin;
2910325fc9f4SBarry Smith   /* currently sctx is memory bleed */
2911325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2912325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2913325fc9f4SBarry Smith   /*
2914325fc9f4SBarry Smith      This should work, but it doesn't
2915325fc9f4SBarry Smith   sctx->ctx = ctx;
2916325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
2917325fc9f4SBarry Smith   */
2918325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
2919cdcf91faSSean Farley   ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
2920325fc9f4SBarry Smith   PetscFunctionReturn(0);
2921325fc9f4SBarry Smith }
2922325fc9f4SBarry Smith 
2923b5b1a830SBarry Smith #undef __FUNCT__
2924b5b1a830SBarry Smith #define __FUNCT__ "TSMonitor_Matlab"
2925b5b1a830SBarry Smith /*
2926b5b1a830SBarry Smith    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
2927b5b1a830SBarry Smith 
2928b5b1a830SBarry Smith    Collective on TS
2929b5b1a830SBarry Smith 
2930b5b1a830SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
2931b5b1a830SBarry Smith @*/
2932cdcf91faSSean Farley PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
2933b5b1a830SBarry Smith {
2934b5b1a830SBarry Smith   PetscErrorCode  ierr;
2935b5b1a830SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
2936a530c242SBarry Smith   int             nlhs = 1,nrhs = 6;
2937b5b1a830SBarry Smith   mxArray         *plhs[1],*prhs[6];
2938b5b1a830SBarry Smith   long long int   lx = 0,ls = 0;
2939b5b1a830SBarry Smith 
2940b5b1a830SBarry Smith   PetscFunctionBegin;
2941cdcf91faSSean Farley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2942b5b1a830SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
2943b5b1a830SBarry Smith 
2944cdcf91faSSean Farley   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
2945b5b1a830SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
2946b5b1a830SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
2947b5b1a830SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)it);
2948b5b1a830SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)time);
2949b5b1a830SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lx);
2950b5b1a830SBarry Smith   prhs[4] =  mxCreateString(sctx->funcname);
2951b5b1a830SBarry Smith   prhs[5] =  sctx->ctx;
2952b5b1a830SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
2953b5b1a830SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
2954b5b1a830SBarry Smith   mxDestroyArray(prhs[0]);
2955b5b1a830SBarry Smith   mxDestroyArray(prhs[1]);
2956b5b1a830SBarry Smith   mxDestroyArray(prhs[2]);
2957b5b1a830SBarry Smith   mxDestroyArray(prhs[3]);
2958b5b1a830SBarry Smith   mxDestroyArray(prhs[4]);
2959b5b1a830SBarry Smith   mxDestroyArray(plhs[0]);
2960b5b1a830SBarry Smith   PetscFunctionReturn(0);
2961b5b1a830SBarry Smith }
2962b5b1a830SBarry Smith 
2963b5b1a830SBarry Smith 
2964b5b1a830SBarry Smith #undef __FUNCT__
2965b5b1a830SBarry Smith #define __FUNCT__ "TSMonitorSetMatlab"
2966b5b1a830SBarry Smith /*
2967b5b1a830SBarry Smith    TSMonitorSetMatlab - Sets the monitor function from Matlab
2968b5b1a830SBarry Smith 
2969b5b1a830SBarry Smith    Level: developer
2970b5b1a830SBarry Smith 
2971b5b1a830SBarry Smith .keywords: TS, nonlinear, set, function
2972b5b1a830SBarry Smith 
2973b5b1a830SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
2974b5b1a830SBarry Smith */
2975cdcf91faSSean Farley PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
2976b5b1a830SBarry Smith {
2977b5b1a830SBarry Smith   PetscErrorCode    ierr;
2978b5b1a830SBarry Smith   TSMatlabContext *sctx;
2979b5b1a830SBarry Smith 
2980b5b1a830SBarry Smith   PetscFunctionBegin;
2981b5b1a830SBarry Smith   /* currently sctx is memory bleed */
2982b5b1a830SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
2983b5b1a830SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
2984b5b1a830SBarry Smith   /*
2985b5b1a830SBarry Smith      This should work, but it doesn't
2986b5b1a830SBarry Smith   sctx->ctx = ctx;
2987b5b1a830SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
2988b5b1a830SBarry Smith   */
2989b5b1a830SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
2990cdcf91faSSean Farley   ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
2991b5b1a830SBarry Smith   PetscFunctionReturn(0);
2992b5b1a830SBarry Smith }
2993325fc9f4SBarry Smith #endif
2994