xref: /petsc/src/ts/interface/ts.c (revision 8caf3d721b8a28a7e6d895d5ad47be9517ebf833)
163dd3a1aSKris Buschelman 
2b45d2f2cSJed Brown #include <petsc-private/tsimpl.h>        /*I "petscts.h"  I*/
3496e6a7aSJed Brown #include <petscdmshell.h>
4d763cef2SBarry Smith 
5d5ba7fb7SMatthew Knepley /* Logging support */
67087cfbeSBarry Smith PetscClassId  TS_CLASSID;
7166c7f25SBarry Smith PetscLogEvent  TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
8d405a339SMatthew Knepley 
94a2ae208SSatish Balay #undef __FUNCT__
10bdad233fSMatthew Knepley #define __FUNCT__ "TSSetTypeFromOptions"
11bdad233fSMatthew Knepley /*
12bdad233fSMatthew Knepley   TSSetTypeFromOptions - Sets the type of ts from user options.
13bdad233fSMatthew Knepley 
14bdad233fSMatthew Knepley   Collective on TS
15bdad233fSMatthew Knepley 
16bdad233fSMatthew Knepley   Input Parameter:
17bdad233fSMatthew Knepley . ts - The ts
18bdad233fSMatthew Knepley 
19bdad233fSMatthew Knepley   Level: intermediate
20bdad233fSMatthew Knepley 
21bdad233fSMatthew Knepley .keywords: TS, set, options, database, type
22bdad233fSMatthew Knepley .seealso: TSSetFromOptions(), TSSetType()
23bdad233fSMatthew Knepley */
246849ba73SBarry Smith static PetscErrorCode TSSetTypeFromOptions(TS ts)
25bdad233fSMatthew Knepley {
26ace3abfcSBarry Smith   PetscBool      opt;
272fc52814SBarry Smith   const char     *defaultType;
28bdad233fSMatthew Knepley   char           typeName[256];
29dfbe8321SBarry Smith   PetscErrorCode ierr;
30bdad233fSMatthew Knepley 
31bdad233fSMatthew Knepley   PetscFunctionBegin;
327adad957SLisandro Dalcin   if (((PetscObject)ts)->type_name) {
337adad957SLisandro Dalcin     defaultType = ((PetscObject)ts)->type_name;
34bdad233fSMatthew Knepley   } else {
359596e0b4SJed Brown     defaultType = TSEULER;
36bdad233fSMatthew Knepley   }
37bdad233fSMatthew Knepley 
38cce0b1b2SLisandro Dalcin   if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
39bdad233fSMatthew Knepley   ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr);
40a7cc72afSBarry Smith   if (opt) {
41bdad233fSMatthew Knepley     ierr = TSSetType(ts, typeName);CHKERRQ(ierr);
42bdad233fSMatthew Knepley   } else {
43bdad233fSMatthew Knepley     ierr = TSSetType(ts, defaultType);CHKERRQ(ierr);
44bdad233fSMatthew Knepley   }
45bdad233fSMatthew Knepley   PetscFunctionReturn(0);
46bdad233fSMatthew Knepley }
47bdad233fSMatthew Knepley 
48bdad233fSMatthew Knepley #undef __FUNCT__
49bdad233fSMatthew Knepley #define __FUNCT__ "TSSetFromOptions"
50bdad233fSMatthew Knepley /*@
51bdad233fSMatthew Knepley    TSSetFromOptions - Sets various TS parameters from user options.
52bdad233fSMatthew Knepley 
53bdad233fSMatthew Knepley    Collective on TS
54bdad233fSMatthew Knepley 
55bdad233fSMatthew Knepley    Input Parameter:
56bdad233fSMatthew Knepley .  ts - the TS context obtained from TSCreate()
57bdad233fSMatthew Knepley 
58bdad233fSMatthew Knepley    Options Database Keys:
594d91e141SJed Brown +  -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP
60bdad233fSMatthew Knepley .  -ts_max_steps maxsteps - maximum number of time-steps to take
613bca7d26SBarry Smith .  -ts_final_time time - maximum time to compute to
62bdad233fSMatthew Knepley .  -ts_dt dt - initial time step
63bdad233fSMatthew Knepley .  -ts_monitor - print information at each timestep
64a6570f20SBarry Smith -  -ts_monitor_draw - plot information at each timestep
65bdad233fSMatthew Knepley 
66bdad233fSMatthew Knepley    Level: beginner
67bdad233fSMatthew Knepley 
68bdad233fSMatthew Knepley .keywords: TS, timestep, set, options, database
69bdad233fSMatthew Knepley 
70a313700dSBarry Smith .seealso: TSGetType()
71bdad233fSMatthew Knepley @*/
727087cfbeSBarry Smith PetscErrorCode  TSSetFromOptions(TS ts)
73bdad233fSMatthew Knepley {
74ace3abfcSBarry Smith   PetscBool      opt,flg;
75dfbe8321SBarry Smith   PetscErrorCode ierr;
76649052a6SBarry Smith   PetscViewer    monviewer;
77eabae89aSBarry Smith   char           monfilename[PETSC_MAX_PATH_LEN];
78089b2837SJed Brown   SNES           snes;
791c3436cfSJed Brown   TSAdapt        adapt;
8031748224SBarry Smith   PetscReal      time_step;
81bdad233fSMatthew Knepley 
82bdad233fSMatthew Knepley   PetscFunctionBegin;
830700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
843194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr);
85a43b19c4SJed Brown     /* Handle TS type options */
86a43b19c4SJed Brown     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
87bdad233fSMatthew Knepley 
88bdad233fSMatthew Knepley     /* Handle generic TS options */
89bdad233fSMatthew Knepley     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
903bca7d26SBarry Smith     ierr = PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
91d8cd7023SBarry Smith     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);CHKERRQ(ierr);
9231748224SBarry Smith     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&time_step,&flg);CHKERRQ(ierr);
9331748224SBarry Smith     if (flg) {
9431748224SBarry Smith       ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
9531748224SBarry Smith     }
96a43b19c4SJed Brown     opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time;
97a43b19c4SJed Brown     ierr = PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);CHKERRQ(ierr);
98461e00b3SSean Farley     if (flg) {ierr = TSSetExactFinalTime(ts,opt);CHKERRQ(ierr);}
99cef5090cSJed Brown     ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr);
100cef5090cSJed Brown     ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr);
101cef5090cSJed Brown     ierr = PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr);
1021c3436cfSJed Brown     ierr = PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,PETSC_NULL);CHKERRQ(ierr);
1031c3436cfSJed Brown     ierr = PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,PETSC_NULL);CHKERRQ(ierr);
104bdad233fSMatthew Knepley 
105bdad233fSMatthew Knepley     /* Monitor options */
106a6570f20SBarry Smith     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
107eabae89aSBarry Smith     if (flg) {
108649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr);
109649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
110bdad233fSMatthew Knepley     }
1115180491cSLisandro Dalcin     ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1125180491cSLisandro Dalcin     if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);}
1135180491cSLisandro Dalcin 
11490d69ab7SBarry Smith     opt  = PETSC_FALSE;
115acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
116a7cc72afSBarry Smith     if (opt) {
117a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
118bdad233fSMatthew Knepley     }
11990d69ab7SBarry Smith     opt  = PETSC_FALSE;
120acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
121a7cc72afSBarry Smith     if (opt) {
1226083293cSBarry Smith       void *ctx;
1236083293cSBarry Smith       ierr = TSMonitorSolutionCreate(ts,PETSC_NULL,&ctx);CHKERRQ(ierr);
1246083293cSBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);CHKERRQ(ierr);
125bdad233fSMatthew Knepley     }
126fb1732b5SBarry Smith     opt  = PETSC_FALSE;
127fb1732b5SBarry Smith     ierr = PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
128fb1732b5SBarry Smith     if (flg) {
129fb1732b5SBarry Smith       PetscViewer ctx;
130fb1732b5SBarry Smith       if (monfilename[0]) {
131fb1732b5SBarry Smith         ierr = PetscViewerBinaryOpen(((PetscObject)ts)->comm,monfilename,FILE_MODE_WRITE,&ctx);CHKERRQ(ierr);
132fb1732b5SBarry Smith       } else {
133fb1732b5SBarry Smith         ctx = PETSC_VIEWER_BINARY_(((PetscObject)ts)->comm);
134fb1732b5SBarry Smith       }
135fb1732b5SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
136fb1732b5SBarry Smith     }
137ed81e22dSJed Brown     opt  = PETSC_FALSE;
138ed81e22dSJed Brown     ierr = PetscOptionsString("-ts_monitor_solution_vtk","Save each time step to a binary file, use filename-%%03D.vts","TSMonitorSolutionVTK",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
139ed81e22dSJed Brown     if (flg) {
140ed81e22dSJed Brown       const char *ptr,*ptr2;
141ed81e22dSJed Brown       char *filetemplate;
142ed81e22dSJed Brown       if (!monfilename[0]) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
143ed81e22dSJed Brown       /* Do some cursory validation of the input. */
144ed81e22dSJed Brown       ierr = PetscStrstr(monfilename,"%",(char**)&ptr);CHKERRQ(ierr);
145ed81e22dSJed Brown       if (!ptr) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
146ed81e22dSJed Brown       for (ptr++ ; ptr && *ptr; ptr++) {
147ed81e22dSJed Brown         ierr = PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);CHKERRQ(ierr);
148ed81e22dSJed Brown         if (!ptr2 && (*ptr < '0' || '9' < *ptr)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03D.vts");
149ed81e22dSJed Brown         if (ptr2) break;
150ed81e22dSJed Brown       }
151ed81e22dSJed Brown       ierr = PetscStrallocpy(monfilename,&filetemplate);CHKERRQ(ierr);
152ed81e22dSJed Brown       ierr = TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);CHKERRQ(ierr);
153ed81e22dSJed Brown     }
154bdad233fSMatthew Knepley 
1551c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
1561c3436cfSJed Brown     ierr = TSAdaptSetFromOptions(adapt);CHKERRQ(ierr);
1571c3436cfSJed Brown 
158d52bd9f3SBarry Smith     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
159d52bd9f3SBarry Smith     if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
160d52bd9f3SBarry Smith 
161bdad233fSMatthew Knepley     /* Handle specific TS options */
162abc0a331SBarry Smith     if (ts->ops->setfromoptions) {
163bdad233fSMatthew Knepley       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
164bdad233fSMatthew Knepley     }
1655d973c19SBarry Smith 
1665d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
1675d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
168bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
169bdad233fSMatthew Knepley   PetscFunctionReturn(0);
170bdad233fSMatthew Knepley }
171bdad233fSMatthew Knepley 
172bdad233fSMatthew Knepley #undef __FUNCT__
173cdcf91faSSean Farley #undef __FUNCT__
1744a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
175a7a1495cSBarry Smith /*@
1768c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
177a7a1495cSBarry Smith       set with TSSetRHSJacobian().
178a7a1495cSBarry Smith 
179a7a1495cSBarry Smith    Collective on TS and Vec
180a7a1495cSBarry Smith 
181a7a1495cSBarry Smith    Input Parameters:
182316643e7SJed Brown +  ts - the TS context
183a7a1495cSBarry Smith .  t - current timestep
184a7a1495cSBarry Smith -  x - input vector
185a7a1495cSBarry Smith 
186a7a1495cSBarry Smith    Output Parameters:
187a7a1495cSBarry Smith +  A - Jacobian matrix
188a7a1495cSBarry Smith .  B - optional preconditioning matrix
189a7a1495cSBarry Smith -  flag - flag indicating matrix structure
190a7a1495cSBarry Smith 
191a7a1495cSBarry Smith    Notes:
192a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
193a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
194a7a1495cSBarry Smith 
19594b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
196a7a1495cSBarry Smith    flag parameter.
197a7a1495cSBarry Smith 
198a7a1495cSBarry Smith    Level: developer
199a7a1495cSBarry Smith 
200a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
201a7a1495cSBarry Smith 
20294b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
203a7a1495cSBarry Smith @*/
2047087cfbeSBarry Smith PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
205a7a1495cSBarry Smith {
206dfbe8321SBarry Smith   PetscErrorCode ierr;
2070e4ef248SJed Brown   PetscInt       Xstate;
20824989b8cSPeter Brune   DM             dm;
20924989b8cSPeter Brune   TSDM           tsdm;
21024989b8cSPeter Brune   TSRHSJacobian  rhsjacobianfunc;
21124989b8cSPeter Brune   void           *ctx;
21224989b8cSPeter Brune   TSIJacobian    ijacobianfunc;
213a7a1495cSBarry Smith 
214a7a1495cSBarry Smith   PetscFunctionBegin;
2150700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2160700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
217c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
21824989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
21924989b8cSPeter Brune   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
22024989b8cSPeter Brune   ierr = DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);CHKERRQ(ierr);
22124989b8cSPeter Brune   ierr = DMTSGetIJacobian(dm,&ijacobianfunc,PETSC_NULL);CHKERRQ(ierr);
2220e4ef248SJed Brown   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
2230e4ef248SJed Brown   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) {
2240e4ef248SJed Brown     *flg = ts->rhsjacobian.mstructure;
2250e4ef248SJed Brown     PetscFunctionReturn(0);
226f8ede8e7SJed Brown   }
227d90be118SSean Farley 
22824989b8cSPeter Brune   if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
229d90be118SSean Farley 
23024989b8cSPeter Brune   if (rhsjacobianfunc) {
231d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
232a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
233a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
23424989b8cSPeter Brune     ierr = (*rhsjacobianfunc)(ts,t,X,A,B,flg,ctx);CHKERRQ(ierr);
235a7a1495cSBarry Smith     PetscStackPop;
236d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
237a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2380700a824SBarry Smith     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
2390700a824SBarry Smith     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
240ef66eb69SBarry Smith   } else {
241214bc6a2SJed Brown     ierr = MatZeroEntries(*A);CHKERRQ(ierr);
242214bc6a2SJed Brown     if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);}
243214bc6a2SJed Brown     *flg = SAME_NONZERO_PATTERN;
244ef66eb69SBarry Smith   }
2450e4ef248SJed Brown   ts->rhsjacobian.time = t;
2460e4ef248SJed Brown   ts->rhsjacobian.X = X;
2470e4ef248SJed Brown   ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr);
2480e4ef248SJed Brown   ts->rhsjacobian.mstructure = *flg;
249a7a1495cSBarry Smith   PetscFunctionReturn(0);
250a7a1495cSBarry Smith }
251a7a1495cSBarry Smith 
2524a2ae208SSatish Balay #undef __FUNCT__
2534a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
254316643e7SJed Brown /*@
255d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
256d763cef2SBarry Smith 
257316643e7SJed Brown    Collective on TS and Vec
258316643e7SJed Brown 
259316643e7SJed Brown    Input Parameters:
260316643e7SJed Brown +  ts - the TS context
261316643e7SJed Brown .  t - current time
262316643e7SJed Brown -  x - state vector
263316643e7SJed Brown 
264316643e7SJed Brown    Output Parameter:
265316643e7SJed Brown .  y - right hand side
266316643e7SJed Brown 
267316643e7SJed Brown    Note:
268316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
269316643e7SJed Brown    is used internally within the nonlinear solvers.
270316643e7SJed Brown 
271316643e7SJed Brown    Level: developer
272316643e7SJed Brown 
273316643e7SJed Brown .keywords: TS, compute
274316643e7SJed Brown 
275316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction()
276316643e7SJed Brown @*/
277dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
278d763cef2SBarry Smith {
279dfbe8321SBarry Smith   PetscErrorCode ierr;
280d763cef2SBarry Smith 
28124989b8cSPeter Brune   TSRHSFunction rhsfunction;
28224989b8cSPeter Brune   TSIFunction   ifunction;
28324989b8cSPeter Brune   void          *ctx;
28424989b8cSPeter Brune   DM            dm;
28524989b8cSPeter Brune 
286d763cef2SBarry Smith   PetscFunctionBegin;
2870700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2880700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2890700a824SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
29024989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
29124989b8cSPeter Brune   ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctx);CHKERRQ(ierr);
29224989b8cSPeter Brune   ierr = DMTSGetIFunction(dm,&ifunction,PETSC_NULL);CHKERRQ(ierr);
293d763cef2SBarry Smith 
29424989b8cSPeter Brune   if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
295d763cef2SBarry Smith 
296d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
29724989b8cSPeter Brune   if (rhsfunction) {
298d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
29924989b8cSPeter Brune     ierr = (*rhsfunction)(ts,t,x,y,ctx);CHKERRQ(ierr);
300d763cef2SBarry Smith     PetscStackPop;
301214bc6a2SJed Brown   } else {
302214bc6a2SJed Brown     ierr = VecZeroEntries(y);CHKERRQ(ierr);
303b2cd27e8SJed Brown   }
30444a41b28SSean Farley 
305d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
306d763cef2SBarry Smith   PetscFunctionReturn(0);
307d763cef2SBarry Smith }
308d763cef2SBarry Smith 
3094a2ae208SSatish Balay #undef __FUNCT__
310214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSVec_Private"
311214bc6a2SJed Brown static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
312214bc6a2SJed Brown {
3132dd45cf8SJed Brown   Vec            F;
314214bc6a2SJed Brown   PetscErrorCode ierr;
315214bc6a2SJed Brown 
316214bc6a2SJed Brown   PetscFunctionBegin;
3170da9a4f0SJed Brown   *Frhs = PETSC_NULL;
3182dd45cf8SJed Brown   ierr = TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
319214bc6a2SJed Brown   if (!ts->Frhs) {
3202dd45cf8SJed Brown     ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr);
321214bc6a2SJed Brown   }
322214bc6a2SJed Brown   *Frhs = ts->Frhs;
323214bc6a2SJed Brown   PetscFunctionReturn(0);
324214bc6a2SJed Brown }
325214bc6a2SJed Brown 
326214bc6a2SJed Brown #undef __FUNCT__
327214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSMats_Private"
328214bc6a2SJed Brown static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
329214bc6a2SJed Brown {
330214bc6a2SJed Brown   Mat            A,B;
3312dd45cf8SJed Brown   PetscErrorCode ierr;
332214bc6a2SJed Brown 
333214bc6a2SJed Brown   PetscFunctionBegin;
334214bc6a2SJed Brown   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
335214bc6a2SJed Brown   if (Arhs) {
336214bc6a2SJed Brown     if (!ts->Arhs) {
337214bc6a2SJed Brown       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
338214bc6a2SJed Brown     }
339214bc6a2SJed Brown     *Arhs = ts->Arhs;
340214bc6a2SJed Brown   }
341214bc6a2SJed Brown   if (Brhs) {
342214bc6a2SJed Brown     if (!ts->Brhs) {
343214bc6a2SJed Brown       ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
344214bc6a2SJed Brown     }
345214bc6a2SJed Brown     *Brhs = ts->Brhs;
346214bc6a2SJed Brown   }
347214bc6a2SJed Brown   PetscFunctionReturn(0);
348214bc6a2SJed Brown }
349214bc6a2SJed Brown 
350214bc6a2SJed Brown #undef __FUNCT__
351316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction"
352316643e7SJed Brown /*@
353316643e7SJed Brown    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
354316643e7SJed Brown 
355316643e7SJed Brown    Collective on TS and Vec
356316643e7SJed Brown 
357316643e7SJed Brown    Input Parameters:
358316643e7SJed Brown +  ts - the TS context
359316643e7SJed Brown .  t - current time
360316643e7SJed Brown .  X - state vector
361214bc6a2SJed Brown .  Xdot - time derivative of state vector
362214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
363316643e7SJed Brown 
364316643e7SJed Brown    Output Parameter:
365316643e7SJed Brown .  Y - right hand side
366316643e7SJed Brown 
367316643e7SJed Brown    Note:
368316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
369316643e7SJed Brown    is used internally within the nonlinear solvers.
370316643e7SJed Brown 
371316643e7SJed Brown    If the user did did not write their equations in implicit form, this
372316643e7SJed Brown    function recasts them in implicit form.
373316643e7SJed Brown 
374316643e7SJed Brown    Level: developer
375316643e7SJed Brown 
376316643e7SJed Brown .keywords: TS, compute
377316643e7SJed Brown 
378316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction()
379316643e7SJed Brown @*/
380214bc6a2SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
381316643e7SJed Brown {
382316643e7SJed Brown   PetscErrorCode ierr;
38324989b8cSPeter Brune   TSIFunction    ifunction;
38424989b8cSPeter Brune   TSRHSFunction  rhsfunction;
38524989b8cSPeter Brune   void           *ctx;
38624989b8cSPeter Brune   DM             dm;
387316643e7SJed Brown 
388316643e7SJed Brown   PetscFunctionBegin;
3890700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3900700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
3910700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
3920700a824SBarry Smith   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
393316643e7SJed Brown 
39424989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
39524989b8cSPeter Brune   ierr = DMTSGetIFunction(dm,&ifunction,&ctx);CHKERRQ(ierr);
39624989b8cSPeter Brune   ierr = DMTSGetRHSFunction(dm,&rhsfunction,PETSC_NULL);CHKERRQ(ierr);
39724989b8cSPeter Brune 
39824989b8cSPeter Brune   if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
399d90be118SSean Farley 
400316643e7SJed Brown   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
40124989b8cSPeter Brune   if (ifunction) {
402316643e7SJed Brown     PetscStackPush("TS user implicit function");
40324989b8cSPeter Brune     ierr = (*ifunction)(ts,t,X,Xdot,Y,ctx);CHKERRQ(ierr);
404316643e7SJed Brown     PetscStackPop;
405214bc6a2SJed Brown   }
406214bc6a2SJed Brown   if (imex) {
40724989b8cSPeter Brune     if (!ifunction) {
4082dd45cf8SJed Brown       ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);
4092dd45cf8SJed Brown     }
41024989b8cSPeter Brune   } else if (rhsfunction) {
41124989b8cSPeter Brune     if (ifunction) {
412214bc6a2SJed Brown       Vec Frhs;
413214bc6a2SJed Brown       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
414214bc6a2SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr);
415214bc6a2SJed Brown       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
4162dd45cf8SJed Brown     } else {
4172dd45cf8SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr);
4182dd45cf8SJed Brown       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
419316643e7SJed Brown     }
4204a6899ffSJed Brown   }
421316643e7SJed Brown   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
422316643e7SJed Brown   PetscFunctionReturn(0);
423316643e7SJed Brown }
424316643e7SJed Brown 
425316643e7SJed Brown #undef __FUNCT__
426316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian"
427316643e7SJed Brown /*@
428316643e7SJed Brown    TSComputeIJacobian - Evaluates the Jacobian of the DAE
429316643e7SJed Brown 
430316643e7SJed Brown    Collective on TS and Vec
431316643e7SJed Brown 
432316643e7SJed Brown    Input
433316643e7SJed Brown       Input Parameters:
434316643e7SJed Brown +  ts - the TS context
435316643e7SJed Brown .  t - current timestep
436316643e7SJed Brown .  X - state vector
437316643e7SJed Brown .  Xdot - time derivative of state vector
438214bc6a2SJed Brown .  shift - shift to apply, see note below
439214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
440316643e7SJed Brown 
441316643e7SJed Brown    Output Parameters:
442316643e7SJed Brown +  A - Jacobian matrix
443316643e7SJed Brown .  B - optional preconditioning matrix
444316643e7SJed Brown -  flag - flag indicating matrix structure
445316643e7SJed Brown 
446316643e7SJed Brown    Notes:
447316643e7SJed Brown    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
448316643e7SJed Brown 
449316643e7SJed Brown    dF/dX + shift*dF/dXdot
450316643e7SJed Brown 
451316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
452316643e7SJed Brown    is used internally within the nonlinear solvers.
453316643e7SJed Brown 
454316643e7SJed Brown    Level: developer
455316643e7SJed Brown 
456316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix
457316643e7SJed Brown 
458316643e7SJed Brown .seealso:  TSSetIJacobian()
45963495f91SJed Brown @*/
460214bc6a2SJed Brown PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
461316643e7SJed Brown {
4620026cea9SSean Farley   PetscInt Xstate, Xdotstate;
463316643e7SJed Brown   PetscErrorCode ierr;
46424989b8cSPeter Brune   TSIJacobian    ijacobian;
46524989b8cSPeter Brune   TSRHSJacobian  rhsjacobian;
46624989b8cSPeter Brune   DM             dm;
46724989b8cSPeter Brune   void           *ctx;
468316643e7SJed Brown 
469316643e7SJed Brown   PetscFunctionBegin;
47024989b8cSPeter Brune 
4710700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4720700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
4730700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
474316643e7SJed Brown   PetscValidPointer(A,6);
4750700a824SBarry Smith   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
476316643e7SJed Brown   PetscValidPointer(B,7);
4770700a824SBarry Smith   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
478316643e7SJed Brown   PetscValidPointer(flg,8);
47924989b8cSPeter Brune 
48024989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
48124989b8cSPeter Brune   ierr = DMTSGetIJacobian(dm,&ijacobian,&ctx);CHKERRQ(ierr);
48224989b8cSPeter Brune   ierr = DMTSGetRHSJacobian(dm,&rhsjacobian,PETSC_NULL);CHKERRQ(ierr);
48324989b8cSPeter Brune 
4840026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
4850026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr);
4860026cea9SSean 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))) {
4870026cea9SSean Farley     *flg = ts->ijacobian.mstructure;
4880026cea9SSean Farley     ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr);
4890026cea9SSean Farley     PetscFunctionReturn(0);
4900026cea9SSean Farley   }
491316643e7SJed Brown 
49224989b8cSPeter Brune   if (!rhsjacobian && !ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
493316643e7SJed Brown 
4944e684422SJed Brown   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
495316643e7SJed Brown   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
49624989b8cSPeter Brune   if (ijacobian) {
4972dd45cf8SJed Brown     *flg = DIFFERENT_NONZERO_PATTERN;
498316643e7SJed Brown     PetscStackPush("TS user implicit Jacobian");
49924989b8cSPeter Brune     ierr = (*ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ctx);CHKERRQ(ierr);
500316643e7SJed Brown     PetscStackPop;
501214bc6a2SJed Brown     /* make sure user returned a correct Jacobian and preconditioner */
502214bc6a2SJed Brown     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
503214bc6a2SJed Brown     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
5044a6899ffSJed Brown   }
505214bc6a2SJed Brown   if (imex) {
50624989b8cSPeter Brune     if (!ijacobian) {  /* system was written as Xdot = F(t,X) */
507214bc6a2SJed Brown       ierr = MatZeroEntries(*A);CHKERRQ(ierr);
5082dd45cf8SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
509214bc6a2SJed Brown       if (*A != *B) {
510214bc6a2SJed Brown         ierr = MatZeroEntries(*B);CHKERRQ(ierr);
511214bc6a2SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
512214bc6a2SJed Brown       }
513214bc6a2SJed Brown       *flg = SAME_PRECONDITIONER;
514214bc6a2SJed Brown     }
515214bc6a2SJed Brown   } else {
51624989b8cSPeter Brune     if (!ijacobian) {
517214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr);
518214bc6a2SJed Brown       ierr = MatScale(*A,-1);CHKERRQ(ierr);
519214bc6a2SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
520316643e7SJed Brown       if (*A != *B) {
521316643e7SJed Brown         ierr = MatScale(*B,-1);CHKERRQ(ierr);
522316643e7SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
523316643e7SJed Brown       }
52424989b8cSPeter Brune     } else if (rhsjacobian) {
525214bc6a2SJed Brown       Mat Arhs,Brhs;
526214bc6a2SJed Brown       MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
527214bc6a2SJed Brown       ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
528214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
529214bc6a2SJed Brown       axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
530214bc6a2SJed Brown       ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr);
531214bc6a2SJed Brown       if (*A != *B) {
532214bc6a2SJed Brown         ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr);
533214bc6a2SJed Brown       }
534214bc6a2SJed Brown       *flg = PetscMin(*flg,flg2);
535214bc6a2SJed Brown     }
536316643e7SJed Brown   }
5370026cea9SSean Farley 
5380026cea9SSean Farley   ts->ijacobian.time = t;
5390026cea9SSean Farley   ts->ijacobian.X = X;
5400026cea9SSean Farley   ts->ijacobian.Xdot = Xdot;
5410026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr);
5420026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr);
5430026cea9SSean Farley   ts->ijacobian.shift = shift;
5440026cea9SSean Farley   ts->ijacobian.imex = imex;
5450026cea9SSean Farley   ts->ijacobian.mstructure = *flg;
546316643e7SJed Brown   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
547316643e7SJed Brown   PetscFunctionReturn(0);
548316643e7SJed Brown }
549316643e7SJed Brown 
550316643e7SJed Brown #undef __FUNCT__
5514a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
552d763cef2SBarry Smith /*@C
553d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
554d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
555d763cef2SBarry Smith 
5563f9fe445SBarry Smith     Logically Collective on TS
557d763cef2SBarry Smith 
558d763cef2SBarry Smith     Input Parameters:
559d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
560ca94891dSJed Brown .   r - vector to put the computed right hand side (or PETSC_NULL to have it created)
561d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
562d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
563d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
564d763cef2SBarry Smith 
565d763cef2SBarry Smith     Calling sequence of func:
56687828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
567d763cef2SBarry Smith 
568d763cef2SBarry Smith +   t - current timestep
569d763cef2SBarry Smith .   u - input vector
570d763cef2SBarry Smith .   F - function vector
571d763cef2SBarry Smith -   ctx - [optional] user-defined function context
572d763cef2SBarry Smith 
573d763cef2SBarry Smith     Level: beginner
574d763cef2SBarry Smith 
575d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
576d763cef2SBarry Smith 
577d6cbdb99SBarry Smith .seealso: TSSetRHSJacobian(), TSSetIJacobian()
578d763cef2SBarry Smith @*/
579089b2837SJed Brown PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
580d763cef2SBarry Smith {
581089b2837SJed Brown   PetscErrorCode ierr;
582089b2837SJed Brown   SNES           snes;
583e856ceecSJed Brown   Vec            ralloc = PETSC_NULL;
58424989b8cSPeter Brune   DM             dm;
585d763cef2SBarry Smith 
586089b2837SJed Brown   PetscFunctionBegin;
5870700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
588ca94891dSJed Brown   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
58924989b8cSPeter Brune 
59024989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
59124989b8cSPeter Brune   ierr = DMTSSetRHSFunction(dm,f,ctx);CHKERRQ(ierr);
592089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
593e856ceecSJed Brown   if (!r && !ts->dm && ts->vec_sol) {
594e856ceecSJed Brown     ierr = VecDuplicate(ts->vec_sol,&ralloc);CHKERRQ(ierr);
595e856ceecSJed Brown     r = ralloc;
596e856ceecSJed Brown   }
597089b2837SJed Brown   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
598e856ceecSJed Brown   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
599d763cef2SBarry Smith   PetscFunctionReturn(0);
600d763cef2SBarry Smith }
601d763cef2SBarry Smith 
6024a2ae208SSatish Balay #undef __FUNCT__
6034a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
604d763cef2SBarry Smith /*@C
605d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
606d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
607d763cef2SBarry Smith 
6083f9fe445SBarry Smith    Logically Collective on TS
609d763cef2SBarry Smith 
610d763cef2SBarry Smith    Input Parameters:
611d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
612d763cef2SBarry Smith .  A   - Jacobian matrix
613d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
614d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
615d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
616d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
617d763cef2SBarry Smith 
618d763cef2SBarry Smith    Calling sequence of func:
61987828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
620d763cef2SBarry Smith 
621d763cef2SBarry Smith +  t - current timestep
622d763cef2SBarry Smith .  u - input vector
623d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
624d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
625d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
62694b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
627d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
628d763cef2SBarry Smith 
629d763cef2SBarry Smith    Notes:
63094b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
631d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
632d763cef2SBarry Smith 
633d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
634d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
63556335db2SHong Zhang    completely new matrix structure (not just different matrix elements)
636d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
637d763cef2SBarry Smith    throughout the global iterations.
638d763cef2SBarry Smith 
639d763cef2SBarry Smith    Level: beginner
640d763cef2SBarry Smith 
641d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
642d763cef2SBarry Smith 
643ab637aeaSJed Brown .seealso: SNESDefaultComputeJacobianColor(), TSSetRHSFunction()
644d763cef2SBarry Smith 
645d763cef2SBarry Smith @*/
646089b2837SJed Brown PetscErrorCode  TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
647d763cef2SBarry Smith {
648277b19d0SLisandro Dalcin   PetscErrorCode ierr;
649089b2837SJed Brown   SNES           snes;
65024989b8cSPeter Brune   DM             dm;
65124989b8cSPeter Brune   TSIJacobian    ijacobian;
652277b19d0SLisandro Dalcin 
653d763cef2SBarry Smith   PetscFunctionBegin;
6540700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
655277b19d0SLisandro Dalcin   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
656277b19d0SLisandro Dalcin   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
657277b19d0SLisandro Dalcin   if (A) PetscCheckSameComm(ts,1,A,2);
658277b19d0SLisandro Dalcin   if (B) PetscCheckSameComm(ts,1,B,3);
659d763cef2SBarry Smith 
66024989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
66124989b8cSPeter Brune   ierr = DMTSSetRHSJacobian(dm,f,ctx);CHKERRQ(ierr);
66224989b8cSPeter Brune   ierr = DMTSGetIJacobian(dm,&ijacobian,PETSC_NULL);CHKERRQ(ierr);
66324989b8cSPeter Brune 
664089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
6655f659677SPeter Brune   if (!ijacobian) {
666089b2837SJed Brown     ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
6670e4ef248SJed Brown   }
6680e4ef248SJed Brown   if (A) {
6690e4ef248SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
6700e4ef248SJed Brown     ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
6710e4ef248SJed Brown     ts->Arhs = A;
6720e4ef248SJed Brown   }
6730e4ef248SJed Brown   if (B) {
6740e4ef248SJed Brown     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
6750e4ef248SJed Brown     ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
6760e4ef248SJed Brown     ts->Brhs = B;
6770e4ef248SJed Brown   }
678d763cef2SBarry Smith   PetscFunctionReturn(0);
679d763cef2SBarry Smith }
680d763cef2SBarry Smith 
681316643e7SJed Brown 
682316643e7SJed Brown #undef __FUNCT__
683316643e7SJed Brown #define __FUNCT__ "TSSetIFunction"
684316643e7SJed Brown /*@C
685316643e7SJed Brown    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
686316643e7SJed Brown 
6873f9fe445SBarry Smith    Logically Collective on TS
688316643e7SJed Brown 
689316643e7SJed Brown    Input Parameters:
690316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
691ca94891dSJed Brown .  r   - vector to hold the residual (or PETSC_NULL to have it created internally)
692316643e7SJed Brown .  f   - the function evaluation routine
693316643e7SJed Brown -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
694316643e7SJed Brown 
695316643e7SJed Brown    Calling sequence of f:
696316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
697316643e7SJed Brown 
698316643e7SJed Brown +  t   - time at step/stage being solved
699316643e7SJed Brown .  u   - state vector
700316643e7SJed Brown .  u_t - time derivative of state vector
701316643e7SJed Brown .  F   - function vector
702316643e7SJed Brown -  ctx - [optional] user-defined context for matrix evaluation routine
703316643e7SJed Brown 
704316643e7SJed Brown    Important:
705d6cbdb99SBarry Smith    The user MUST call either this routine, TSSetRHSFunction().  This routine must be used when not solving an ODE, for example a DAE.
706316643e7SJed Brown 
707316643e7SJed Brown    Level: beginner
708316643e7SJed Brown 
709316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian
710316643e7SJed Brown 
711d6cbdb99SBarry Smith .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
712316643e7SJed Brown @*/
713089b2837SJed Brown PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
714316643e7SJed Brown {
715089b2837SJed Brown   PetscErrorCode ierr;
716089b2837SJed Brown   SNES           snes;
717e856ceecSJed Brown   Vec            resalloc = PETSC_NULL;
71824989b8cSPeter Brune   DM             dm;
719316643e7SJed Brown 
720316643e7SJed Brown   PetscFunctionBegin;
7210700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
722ca94891dSJed Brown   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
72324989b8cSPeter Brune 
72424989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
72524989b8cSPeter Brune   ierr = DMTSSetIFunction(dm,f,ctx);CHKERRQ(ierr);
72624989b8cSPeter Brune 
727089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
728e856ceecSJed Brown   if (!res && !ts->dm && ts->vec_sol) {
729e856ceecSJed Brown     ierr = VecDuplicate(ts->vec_sol,&resalloc);CHKERRQ(ierr);
730e856ceecSJed Brown     res = resalloc;
731e856ceecSJed Brown   }
732089b2837SJed Brown   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
733e856ceecSJed Brown   ierr = VecDestroy(&resalloc);CHKERRQ(ierr);
73424989b8cSPeter Brune 
735089b2837SJed Brown   PetscFunctionReturn(0);
736089b2837SJed Brown }
737089b2837SJed Brown 
738089b2837SJed Brown #undef __FUNCT__
739089b2837SJed Brown #define __FUNCT__ "TSGetIFunction"
740089b2837SJed Brown /*@C
741089b2837SJed Brown    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
742089b2837SJed Brown 
743089b2837SJed Brown    Not Collective
744089b2837SJed Brown 
745089b2837SJed Brown    Input Parameter:
746089b2837SJed Brown .  ts - the TS context
747089b2837SJed Brown 
748089b2837SJed Brown    Output Parameter:
749089b2837SJed Brown +  r - vector to hold residual (or PETSC_NULL)
750089b2837SJed Brown .  func - the function to compute residual (or PETSC_NULL)
751089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
752089b2837SJed Brown 
753089b2837SJed Brown    Level: advanced
754089b2837SJed Brown 
755089b2837SJed Brown .keywords: TS, nonlinear, get, function
756089b2837SJed Brown 
757089b2837SJed Brown .seealso: TSSetIFunction(), SNESGetFunction()
758089b2837SJed Brown @*/
759089b2837SJed Brown PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
760089b2837SJed Brown {
761089b2837SJed Brown   PetscErrorCode ierr;
762089b2837SJed Brown   SNES snes;
76324989b8cSPeter Brune   DM   dm;
764089b2837SJed Brown 
765089b2837SJed Brown   PetscFunctionBegin;
766089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
767089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
768089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
76924989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
77024989b8cSPeter Brune   ierr = DMTSGetIFunction(dm,func,ctx);CHKERRQ(ierr);
771089b2837SJed Brown   PetscFunctionReturn(0);
772089b2837SJed Brown }
773089b2837SJed Brown 
774089b2837SJed Brown #undef __FUNCT__
775089b2837SJed Brown #define __FUNCT__ "TSGetRHSFunction"
776089b2837SJed Brown /*@C
777089b2837SJed Brown    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
778089b2837SJed Brown 
779089b2837SJed Brown    Not Collective
780089b2837SJed Brown 
781089b2837SJed Brown    Input Parameter:
782089b2837SJed Brown .  ts - the TS context
783089b2837SJed Brown 
784089b2837SJed Brown    Output Parameter:
785089b2837SJed Brown +  r - vector to hold computed right hand side (or PETSC_NULL)
786089b2837SJed Brown .  func - the function to compute right hand side (or PETSC_NULL)
787089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
788089b2837SJed Brown 
789089b2837SJed Brown    Level: advanced
790089b2837SJed Brown 
791089b2837SJed Brown .keywords: TS, nonlinear, get, function
792089b2837SJed Brown 
793089b2837SJed Brown .seealso: TSSetRhsfunction(), SNESGetFunction()
794089b2837SJed Brown @*/
795089b2837SJed Brown PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
796089b2837SJed Brown {
797089b2837SJed Brown   PetscErrorCode ierr;
798089b2837SJed Brown   SNES snes;
79924989b8cSPeter Brune   DM   dm;
800089b2837SJed Brown 
801089b2837SJed Brown   PetscFunctionBegin;
802089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
803089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
804089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
80524989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
80624989b8cSPeter Brune   ierr = DMTSGetRHSFunction(dm,func,ctx);CHKERRQ(ierr);
807316643e7SJed Brown   PetscFunctionReturn(0);
808316643e7SJed Brown }
809316643e7SJed Brown 
810316643e7SJed Brown #undef __FUNCT__
811316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian"
812316643e7SJed Brown /*@C
813a4f0a591SBarry Smith    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
814a4f0a591SBarry Smith         you provided with TSSetIFunction().
815316643e7SJed Brown 
8163f9fe445SBarry Smith    Logically Collective on TS
817316643e7SJed Brown 
818316643e7SJed Brown    Input Parameters:
819316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
820316643e7SJed Brown .  A   - Jacobian matrix
821316643e7SJed Brown .  B   - preconditioning matrix for A (may be same as A)
822316643e7SJed Brown .  f   - the Jacobian evaluation routine
823316643e7SJed Brown -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
824316643e7SJed Brown 
825316643e7SJed Brown    Calling sequence of f:
8261b4a444bSJed Brown $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
827316643e7SJed Brown 
828316643e7SJed Brown +  t    - time at step/stage being solved
8291b4a444bSJed Brown .  U    - state vector
8301b4a444bSJed Brown .  U_t  - time derivative of state vector
831316643e7SJed Brown .  a    - shift
8321b4a444bSJed Brown .  A    - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
833316643e7SJed Brown .  B    - preconditioning matrix for A, may be same as A
834316643e7SJed Brown .  flag - flag indicating information about the preconditioner matrix
835316643e7SJed Brown           structure (same as flag in KSPSetOperators())
836316643e7SJed Brown -  ctx  - [optional] user-defined context for matrix evaluation routine
837316643e7SJed Brown 
838316643e7SJed Brown    Notes:
839316643e7SJed Brown    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
840316643e7SJed Brown 
841a4f0a591SBarry Smith    The matrix dF/dU + a*dF/dU_t you provide turns out to be
842a4f0a591SBarry 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.
843a4f0a591SBarry Smith    The time integrator internally approximates U_t by W+a*U where the positive "shift"
844a4f0a591SBarry Smith    a and vector W depend on the integration method, step size, and past states. For example with
845a4f0a591SBarry Smith    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
846a4f0a591SBarry Smith    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
847a4f0a591SBarry Smith 
848316643e7SJed Brown    Level: beginner
849316643e7SJed Brown 
850316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian
851316643e7SJed Brown 
852ab637aeaSJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESDefaultComputeJacobianColor(), SNESDefaultComputeJacobian()
853316643e7SJed Brown 
854316643e7SJed Brown @*/
8557087cfbeSBarry Smith PetscErrorCode  TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
856316643e7SJed Brown {
857316643e7SJed Brown   PetscErrorCode ierr;
858089b2837SJed Brown   SNES           snes;
85924989b8cSPeter Brune   DM             dm;
860316643e7SJed Brown 
861316643e7SJed Brown   PetscFunctionBegin;
8620700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8630700a824SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
8640700a824SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
865316643e7SJed Brown   if (A) PetscCheckSameComm(ts,1,A,2);
866316643e7SJed Brown   if (B) PetscCheckSameComm(ts,1,B,3);
86724989b8cSPeter Brune 
86824989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
86924989b8cSPeter Brune   ierr = DMTSSetIJacobian(dm,f,ctx);CHKERRQ(ierr);
87024989b8cSPeter Brune 
871089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
872089b2837SJed Brown   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
873316643e7SJed Brown   PetscFunctionReturn(0);
874316643e7SJed Brown }
875316643e7SJed Brown 
8764a2ae208SSatish Balay #undef __FUNCT__
8774a2ae208SSatish Balay #define __FUNCT__ "TSView"
8787e2c5f70SBarry Smith /*@C
879d763cef2SBarry Smith     TSView - Prints the TS data structure.
880d763cef2SBarry Smith 
8814c49b128SBarry Smith     Collective on TS
882d763cef2SBarry Smith 
883d763cef2SBarry Smith     Input Parameters:
884d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
885d763cef2SBarry Smith -   viewer - visualization context
886d763cef2SBarry Smith 
887d763cef2SBarry Smith     Options Database Key:
888d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
889d763cef2SBarry Smith 
890d763cef2SBarry Smith     Notes:
891d763cef2SBarry Smith     The available visualization contexts include
892b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
893b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
894d763cef2SBarry Smith          output where only the first processor opens
895d763cef2SBarry Smith          the file.  All other processors send their
896d763cef2SBarry Smith          data to the first processor to print.
897d763cef2SBarry Smith 
898d763cef2SBarry Smith     The user can open an alternative visualization context with
899b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
900d763cef2SBarry Smith 
901d763cef2SBarry Smith     Level: beginner
902d763cef2SBarry Smith 
903d763cef2SBarry Smith .keywords: TS, timestep, view
904d763cef2SBarry Smith 
905b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
906d763cef2SBarry Smith @*/
9077087cfbeSBarry Smith PetscErrorCode  TSView(TS ts,PetscViewer viewer)
908d763cef2SBarry Smith {
909dfbe8321SBarry Smith   PetscErrorCode ierr;
910a313700dSBarry Smith   const TSType   type;
911089b2837SJed Brown   PetscBool      iascii,isstring,isundials;
912d763cef2SBarry Smith 
913d763cef2SBarry Smith   PetscFunctionBegin;
9140700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
9153050cee2SBarry Smith   if (!viewer) {
9167adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
9173050cee2SBarry Smith   }
9180700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
919c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
920fd16b177SBarry Smith 
921251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
922251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
92332077d6dSBarry Smith   if (iascii) {
924317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr);
92577431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
926a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
927d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
9285ef26d82SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);CHKERRQ(ierr);
929c610991cSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr);
930d763cef2SBarry Smith     }
9315ef26d82SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);CHKERRQ(ierr);
932193ac0bcSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
933d52bd9f3SBarry Smith     if (ts->ops->view) {
934d52bd9f3SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
935d52bd9f3SBarry Smith       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
936d52bd9f3SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
937d52bd9f3SBarry Smith     }
9380f5bd95cSBarry Smith   } else if (isstring) {
939a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
940b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
941d763cef2SBarry Smith   }
942b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
943251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
944b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
945d763cef2SBarry Smith   PetscFunctionReturn(0);
946d763cef2SBarry Smith }
947d763cef2SBarry Smith 
948d763cef2SBarry Smith 
9494a2ae208SSatish Balay #undef __FUNCT__
9504a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
951b07ff414SBarry Smith /*@
952d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
953d763cef2SBarry Smith    the timesteppers.
954d763cef2SBarry Smith 
9553f9fe445SBarry Smith    Logically Collective on TS
956d763cef2SBarry Smith 
957d763cef2SBarry Smith    Input Parameters:
958d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
959d763cef2SBarry Smith -  usrP - optional user context
960d763cef2SBarry Smith 
961d763cef2SBarry Smith    Level: intermediate
962d763cef2SBarry Smith 
963d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
964d763cef2SBarry Smith 
965d763cef2SBarry Smith .seealso: TSGetApplicationContext()
966d763cef2SBarry Smith @*/
9677087cfbeSBarry Smith PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
968d763cef2SBarry Smith {
969d763cef2SBarry Smith   PetscFunctionBegin;
9700700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
971d763cef2SBarry Smith   ts->user = usrP;
972d763cef2SBarry Smith   PetscFunctionReturn(0);
973d763cef2SBarry Smith }
974d763cef2SBarry Smith 
9754a2ae208SSatish Balay #undef __FUNCT__
9764a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
977b07ff414SBarry Smith /*@
978d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
979d763cef2SBarry Smith     timestepper.
980d763cef2SBarry Smith 
981d763cef2SBarry Smith     Not Collective
982d763cef2SBarry Smith 
983d763cef2SBarry Smith     Input Parameter:
984d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
985d763cef2SBarry Smith 
986d763cef2SBarry Smith     Output Parameter:
987d763cef2SBarry Smith .   usrP - user context
988d763cef2SBarry Smith 
989d763cef2SBarry Smith     Level: intermediate
990d763cef2SBarry Smith 
991d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
992d763cef2SBarry Smith 
993d763cef2SBarry Smith .seealso: TSSetApplicationContext()
994d763cef2SBarry Smith @*/
995e71120c6SJed Brown PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
996d763cef2SBarry Smith {
997d763cef2SBarry Smith   PetscFunctionBegin;
9980700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
999e71120c6SJed Brown   *(void**)usrP = ts->user;
1000d763cef2SBarry Smith   PetscFunctionReturn(0);
1001d763cef2SBarry Smith }
1002d763cef2SBarry Smith 
10034a2ae208SSatish Balay #undef __FUNCT__
10044a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
1005d763cef2SBarry Smith /*@
1006b8123daeSJed Brown    TSGetTimeStepNumber - Gets the number of time steps completed.
1007d763cef2SBarry Smith 
1008d763cef2SBarry Smith    Not Collective
1009d763cef2SBarry Smith 
1010d763cef2SBarry Smith    Input Parameter:
1011d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1012d763cef2SBarry Smith 
1013d763cef2SBarry Smith    Output Parameter:
1014b8123daeSJed Brown .  iter - number of steps completed so far
1015d763cef2SBarry Smith 
1016d763cef2SBarry Smith    Level: intermediate
1017d763cef2SBarry Smith 
1018d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
1019b8123daeSJed Brown .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStep()
1020d763cef2SBarry Smith @*/
10217087cfbeSBarry Smith PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
1022d763cef2SBarry Smith {
1023d763cef2SBarry Smith   PetscFunctionBegin;
10240700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10254482741eSBarry Smith   PetscValidIntPointer(iter,2);
1026d763cef2SBarry Smith   *iter = ts->steps;
1027d763cef2SBarry Smith   PetscFunctionReturn(0);
1028d763cef2SBarry Smith }
1029d763cef2SBarry Smith 
10304a2ae208SSatish Balay #undef __FUNCT__
10314a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
1032d763cef2SBarry Smith /*@
1033d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
1034d763cef2SBarry Smith    as well as the initial time.
1035d763cef2SBarry Smith 
10363f9fe445SBarry Smith    Logically Collective on TS
1037d763cef2SBarry Smith 
1038d763cef2SBarry Smith    Input Parameters:
1039d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1040d763cef2SBarry Smith .  initial_time - the initial time
1041d763cef2SBarry Smith -  time_step - the size of the timestep
1042d763cef2SBarry Smith 
1043d763cef2SBarry Smith    Level: intermediate
1044d763cef2SBarry Smith 
1045d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
1046d763cef2SBarry Smith 
1047d763cef2SBarry Smith .keywords: TS, set, initial, timestep
1048d763cef2SBarry Smith @*/
10497087cfbeSBarry Smith PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
1050d763cef2SBarry Smith {
1051e144a568SJed Brown   PetscErrorCode ierr;
1052e144a568SJed Brown 
1053d763cef2SBarry Smith   PetscFunctionBegin;
10540700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1055e144a568SJed Brown   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
1056d8cd7023SBarry Smith   ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr);
1057d763cef2SBarry Smith   PetscFunctionReturn(0);
1058d763cef2SBarry Smith }
1059d763cef2SBarry Smith 
10604a2ae208SSatish Balay #undef __FUNCT__
10614a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
1062d763cef2SBarry Smith /*@
1063d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
1064d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
1065d763cef2SBarry Smith 
10663f9fe445SBarry Smith    Logically Collective on TS
1067d763cef2SBarry Smith 
1068d763cef2SBarry Smith    Input Parameters:
1069d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1070d763cef2SBarry Smith -  time_step - the size of the timestep
1071d763cef2SBarry Smith 
1072d763cef2SBarry Smith    Level: intermediate
1073d763cef2SBarry Smith 
1074d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1075d763cef2SBarry Smith 
1076d763cef2SBarry Smith .keywords: TS, set, timestep
1077d763cef2SBarry Smith @*/
10787087cfbeSBarry Smith PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
1079d763cef2SBarry Smith {
1080d763cef2SBarry Smith   PetscFunctionBegin;
10810700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1082c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,time_step,2);
1083d763cef2SBarry Smith   ts->time_step = time_step;
108431748224SBarry Smith   ts->time_step_orig = time_step;
1085d763cef2SBarry Smith   PetscFunctionReturn(0);
1086d763cef2SBarry Smith }
1087d763cef2SBarry Smith 
10884a2ae208SSatish Balay #undef __FUNCT__
1089a43b19c4SJed Brown #define __FUNCT__ "TSSetExactFinalTime"
1090a43b19c4SJed Brown /*@
1091a43b19c4SJed Brown    TSSetExactFinalTime - Determines whether to interpolate solution to the
1092a43b19c4SJed Brown       exact final time requested by the user or just returns it at the final time
1093a43b19c4SJed Brown       it computed.
1094a43b19c4SJed Brown 
1095a43b19c4SJed Brown   Logically Collective on TS
1096a43b19c4SJed Brown 
1097a43b19c4SJed Brown    Input Parameter:
1098a43b19c4SJed Brown +   ts - the time-step context
1099a43b19c4SJed Brown -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
1100a43b19c4SJed Brown 
1101a43b19c4SJed Brown    Level: beginner
1102a43b19c4SJed Brown 
1103a43b19c4SJed Brown .seealso: TSSetDuration()
1104a43b19c4SJed Brown @*/
1105a43b19c4SJed Brown PetscErrorCode  TSSetExactFinalTime(TS ts,PetscBool flg)
1106a43b19c4SJed Brown {
1107a43b19c4SJed Brown 
1108a43b19c4SJed Brown   PetscFunctionBegin;
1109a43b19c4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1110d8cd7023SBarry Smith   PetscValidLogicalCollectiveBool(ts,flg,2);
1111a43b19c4SJed Brown   ts->exact_final_time = flg;
1112a43b19c4SJed Brown   PetscFunctionReturn(0);
1113a43b19c4SJed Brown }
1114a43b19c4SJed Brown 
1115a43b19c4SJed Brown #undef __FUNCT__
11164a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
1117d763cef2SBarry Smith /*@
1118d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
1119d763cef2SBarry Smith 
1120d763cef2SBarry Smith    Not Collective
1121d763cef2SBarry Smith 
1122d763cef2SBarry Smith    Input Parameter:
1123d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1124d763cef2SBarry Smith 
1125d763cef2SBarry Smith    Output Parameter:
1126d763cef2SBarry Smith .  dt - the current timestep size
1127d763cef2SBarry Smith 
1128d763cef2SBarry Smith    Level: intermediate
1129d763cef2SBarry Smith 
1130d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1131d763cef2SBarry Smith 
1132d763cef2SBarry Smith .keywords: TS, get, timestep
1133d763cef2SBarry Smith @*/
11347087cfbeSBarry Smith PetscErrorCode  TSGetTimeStep(TS ts,PetscReal* dt)
1135d763cef2SBarry Smith {
1136d763cef2SBarry Smith   PetscFunctionBegin;
11370700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1138f7cf8827SBarry Smith   PetscValidRealPointer(dt,2);
1139d763cef2SBarry Smith   *dt = ts->time_step;
1140d763cef2SBarry Smith   PetscFunctionReturn(0);
1141d763cef2SBarry Smith }
1142d763cef2SBarry Smith 
11434a2ae208SSatish Balay #undef __FUNCT__
11444a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
1145d8e5e3e6SSatish Balay /*@
1146d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
1147d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
1148d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
1149d763cef2SBarry Smith    the solution at the next timestep has been calculated.
1150d763cef2SBarry Smith 
1151d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
1152d763cef2SBarry Smith 
1153d763cef2SBarry Smith    Input Parameter:
1154d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1155d763cef2SBarry Smith 
1156d763cef2SBarry Smith    Output Parameter:
1157d763cef2SBarry Smith .  v - the vector containing the solution
1158d763cef2SBarry Smith 
1159d763cef2SBarry Smith    Level: intermediate
1160d763cef2SBarry Smith 
1161d763cef2SBarry Smith .seealso: TSGetTimeStep()
1162d763cef2SBarry Smith 
1163d763cef2SBarry Smith .keywords: TS, timestep, get, solution
1164d763cef2SBarry Smith @*/
11657087cfbeSBarry Smith PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1166d763cef2SBarry Smith {
1167d763cef2SBarry Smith   PetscFunctionBegin;
11680700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
11694482741eSBarry Smith   PetscValidPointer(v,2);
11708737fe31SLisandro Dalcin   *v = ts->vec_sol;
1171d763cef2SBarry Smith   PetscFunctionReturn(0);
1172d763cef2SBarry Smith }
1173d763cef2SBarry Smith 
1174bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
11754a2ae208SSatish Balay #undef __FUNCT__
1176bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
1177d8e5e3e6SSatish Balay /*@
1178bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
1179d763cef2SBarry Smith 
1180bdad233fSMatthew Knepley   Not collective
1181d763cef2SBarry Smith 
1182bdad233fSMatthew Knepley   Input Parameters:
1183bdad233fSMatthew Knepley + ts   - The TS
1184bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1185d763cef2SBarry Smith .vb
1186d763cef2SBarry Smith          U_t = A U
1187d763cef2SBarry Smith          U_t = A(t) U
1188d763cef2SBarry Smith          U_t = F(t,U)
1189d763cef2SBarry Smith .ve
1190d763cef2SBarry Smith 
1191d763cef2SBarry Smith    Level: beginner
1192d763cef2SBarry Smith 
1193bdad233fSMatthew Knepley .keywords: TS, problem type
1194bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1195d763cef2SBarry Smith @*/
11967087cfbeSBarry Smith PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1197a7cc72afSBarry Smith {
11989e2a6581SJed Brown   PetscErrorCode ierr;
11999e2a6581SJed Brown 
1200d763cef2SBarry Smith   PetscFunctionBegin;
12010700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1202bdad233fSMatthew Knepley   ts->problem_type = type;
12039e2a6581SJed Brown   if (type == TS_LINEAR) {
12049e2a6581SJed Brown     SNES snes;
12059e2a6581SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
12069e2a6581SJed Brown     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
12079e2a6581SJed Brown   }
1208d763cef2SBarry Smith   PetscFunctionReturn(0);
1209d763cef2SBarry Smith }
1210d763cef2SBarry Smith 
1211bdad233fSMatthew Knepley #undef __FUNCT__
1212bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
1213bdad233fSMatthew Knepley /*@C
1214bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
1215bdad233fSMatthew Knepley 
1216bdad233fSMatthew Knepley   Not collective
1217bdad233fSMatthew Knepley 
1218bdad233fSMatthew Knepley   Input Parameter:
1219bdad233fSMatthew Knepley . ts   - The TS
1220bdad233fSMatthew Knepley 
1221bdad233fSMatthew Knepley   Output Parameter:
1222bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1223bdad233fSMatthew Knepley .vb
1224089b2837SJed Brown          M U_t = A U
1225089b2837SJed Brown          M(t) U_t = A(t) U
1226bdad233fSMatthew Knepley          U_t = F(t,U)
1227bdad233fSMatthew Knepley .ve
1228bdad233fSMatthew Knepley 
1229bdad233fSMatthew Knepley    Level: beginner
1230bdad233fSMatthew Knepley 
1231bdad233fSMatthew Knepley .keywords: TS, problem type
1232bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1233bdad233fSMatthew Knepley @*/
12347087cfbeSBarry Smith PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1235a7cc72afSBarry Smith {
1236bdad233fSMatthew Knepley   PetscFunctionBegin;
12370700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
12384482741eSBarry Smith   PetscValidIntPointer(type,2);
1239bdad233fSMatthew Knepley   *type = ts->problem_type;
1240bdad233fSMatthew Knepley   PetscFunctionReturn(0);
1241bdad233fSMatthew Knepley }
1242d763cef2SBarry Smith 
12434a2ae208SSatish Balay #undef __FUNCT__
12444a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
1245d763cef2SBarry Smith /*@
1246d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
1247d763cef2SBarry Smith    of a timestepper.
1248d763cef2SBarry Smith 
1249d763cef2SBarry Smith    Collective on TS
1250d763cef2SBarry Smith 
1251d763cef2SBarry Smith    Input Parameter:
1252d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1253d763cef2SBarry Smith 
1254d763cef2SBarry Smith    Notes:
1255d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
1256d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
1257d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
1258d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
1259d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
1260d763cef2SBarry Smith 
1261d763cef2SBarry Smith    Level: advanced
1262d763cef2SBarry Smith 
1263d763cef2SBarry Smith .keywords: TS, timestep, setup
1264d763cef2SBarry Smith 
1265d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
1266d763cef2SBarry Smith @*/
12677087cfbeSBarry Smith PetscErrorCode  TSSetUp(TS ts)
1268d763cef2SBarry Smith {
1269dfbe8321SBarry Smith   PetscErrorCode ierr;
12706c6b9e74SPeter Brune   DM             dm;
12716c6b9e74SPeter Brune   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
12726c6b9e74SPeter Brune   PetscErrorCode (*jac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
12736c6b9e74SPeter Brune   TSIJacobian    ijac;
12746c6b9e74SPeter Brune   TSRHSJacobian  rhsjac;
1275d763cef2SBarry Smith 
1276d763cef2SBarry Smith   PetscFunctionBegin;
12770700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1278277b19d0SLisandro Dalcin   if (ts->setupcalled) PetscFunctionReturn(0);
1279277b19d0SLisandro Dalcin 
12807adad957SLisandro Dalcin   if (!((PetscObject)ts)->type_name) {
12819596e0b4SJed Brown     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1282d763cef2SBarry Smith   }
1283a43b19c4SJed Brown   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1284277b19d0SLisandro Dalcin 
1285277b19d0SLisandro Dalcin   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1286277b19d0SLisandro Dalcin 
12871c3436cfSJed Brown   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
12881c3436cfSJed Brown 
1289277b19d0SLisandro Dalcin   if (ts->ops->setup) {
1290000e7ae3SMatthew Knepley     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1291277b19d0SLisandro Dalcin   }
1292277b19d0SLisandro Dalcin 
12936c6b9e74SPeter Brune   /* in the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
12946c6b9e74SPeter Brune    to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
12956c6b9e74SPeter Brune    */
12966c6b9e74SPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
12976c6b9e74SPeter Brune   ierr = DMSNESGetFunction(dm,&func,PETSC_NULL);CHKERRQ(ierr);
12986c6b9e74SPeter Brune   if (!func) {
12996c6b9e74SPeter Brune     ierr =DMSNESSetFunction(dm,SNESTSFormFunction,ts);CHKERRQ(ierr);
13006c6b9e74SPeter Brune   }
13016c6b9e74SPeter Brune   /* if the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it.
13026c6b9e74SPeter Brune      Otherwise, the SNES will use coloring internally to form the Jacobian.
13036c6b9e74SPeter Brune    */
13046c6b9e74SPeter Brune   ierr = DMSNESGetJacobian(dm,&jac,PETSC_NULL);CHKERRQ(ierr);
13056c6b9e74SPeter Brune   ierr = DMTSGetIJacobian(dm,&ijac,PETSC_NULL);CHKERRQ(ierr);
13066c6b9e74SPeter Brune   ierr = DMTSGetRHSJacobian(dm,&rhsjac,PETSC_NULL);CHKERRQ(ierr);
13076c6b9e74SPeter Brune   if (!jac && (ijac || rhsjac)) {
13086c6b9e74SPeter Brune     ierr = DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);CHKERRQ(ierr);
13096c6b9e74SPeter Brune   }
1310277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_TRUE;
1311277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1312277b19d0SLisandro Dalcin }
1313277b19d0SLisandro Dalcin 
1314277b19d0SLisandro Dalcin #undef __FUNCT__
1315277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset"
1316277b19d0SLisandro Dalcin /*@
1317277b19d0SLisandro Dalcin    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1318277b19d0SLisandro Dalcin 
1319277b19d0SLisandro Dalcin    Collective on TS
1320277b19d0SLisandro Dalcin 
1321277b19d0SLisandro Dalcin    Input Parameter:
1322277b19d0SLisandro Dalcin .  ts - the TS context obtained from TSCreate()
1323277b19d0SLisandro Dalcin 
1324277b19d0SLisandro Dalcin    Level: beginner
1325277b19d0SLisandro Dalcin 
1326277b19d0SLisandro Dalcin .keywords: TS, timestep, reset
1327277b19d0SLisandro Dalcin 
1328277b19d0SLisandro Dalcin .seealso: TSCreate(), TSSetup(), TSDestroy()
1329277b19d0SLisandro Dalcin @*/
1330277b19d0SLisandro Dalcin PetscErrorCode  TSReset(TS ts)
1331277b19d0SLisandro Dalcin {
1332277b19d0SLisandro Dalcin   PetscErrorCode ierr;
1333277b19d0SLisandro Dalcin 
1334277b19d0SLisandro Dalcin   PetscFunctionBegin;
1335277b19d0SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1336277b19d0SLisandro Dalcin   if (ts->ops->reset) {
1337277b19d0SLisandro Dalcin     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1338277b19d0SLisandro Dalcin   }
1339277b19d0SLisandro Dalcin   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
13404e684422SJed Brown   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
13414e684422SJed Brown   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1342214bc6a2SJed Brown   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
13436bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1344e3d84a46SJed Brown   ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
1345e3d84a46SJed Brown   ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
134638637c2eSJed Brown   ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);
1347277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_FALSE;
1348d763cef2SBarry Smith   PetscFunctionReturn(0);
1349d763cef2SBarry Smith }
1350d763cef2SBarry Smith 
13514a2ae208SSatish Balay #undef __FUNCT__
13524a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
1353d8e5e3e6SSatish Balay /*@
1354d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
1355d763cef2SBarry Smith    with TSCreate().
1356d763cef2SBarry Smith 
1357d763cef2SBarry Smith    Collective on TS
1358d763cef2SBarry Smith 
1359d763cef2SBarry Smith    Input Parameter:
1360d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1361d763cef2SBarry Smith 
1362d763cef2SBarry Smith    Level: beginner
1363d763cef2SBarry Smith 
1364d763cef2SBarry Smith .keywords: TS, timestepper, destroy
1365d763cef2SBarry Smith 
1366d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
1367d763cef2SBarry Smith @*/
13686bf464f9SBarry Smith PetscErrorCode  TSDestroy(TS *ts)
1369d763cef2SBarry Smith {
13706849ba73SBarry Smith   PetscErrorCode ierr;
1371d763cef2SBarry Smith 
1372d763cef2SBarry Smith   PetscFunctionBegin;
13736bf464f9SBarry Smith   if (!*ts) PetscFunctionReturn(0);
13746bf464f9SBarry Smith   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
13756bf464f9SBarry Smith   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1376d763cef2SBarry Smith 
13776bf464f9SBarry Smith   ierr = TSReset((*ts));CHKERRQ(ierr);
1378277b19d0SLisandro Dalcin 
1379be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
13806bf464f9SBarry Smith   ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr);
13816bf464f9SBarry Smith   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
13826d4c513bSLisandro Dalcin 
138384df9cb4SJed Brown   ierr = TSAdaptDestroy(&(*ts)->adapt);CHKERRQ(ierr);
13846bf464f9SBarry Smith   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
13856bf464f9SBarry Smith   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
13866bf464f9SBarry Smith   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
13876d4c513bSLisandro Dalcin 
1388a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1389d763cef2SBarry Smith   PetscFunctionReturn(0);
1390d763cef2SBarry Smith }
1391d763cef2SBarry Smith 
13924a2ae208SSatish Balay #undef __FUNCT__
13934a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
1394d8e5e3e6SSatish Balay /*@
1395d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1396d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
1397d763cef2SBarry Smith 
1398d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
1399d763cef2SBarry Smith 
1400d763cef2SBarry Smith    Input Parameter:
1401d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1402d763cef2SBarry Smith 
1403d763cef2SBarry Smith    Output Parameter:
1404d763cef2SBarry Smith .  snes - the nonlinear solver context
1405d763cef2SBarry Smith 
1406d763cef2SBarry Smith    Notes:
1407d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
1408d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
140994b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
1410d763cef2SBarry Smith 
1411d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
1412d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
1413d763cef2SBarry Smith 
1414d763cef2SBarry Smith    Level: beginner
1415d763cef2SBarry Smith 
1416d763cef2SBarry Smith .keywords: timestep, get, SNES
1417d763cef2SBarry Smith @*/
14187087cfbeSBarry Smith PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1419d763cef2SBarry Smith {
1420d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1421d372ba47SLisandro Dalcin 
1422d763cef2SBarry Smith   PetscFunctionBegin;
14230700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
14244482741eSBarry Smith   PetscValidPointer(snes,2);
1425d372ba47SLisandro Dalcin   if (!ts->snes) {
1426d372ba47SLisandro Dalcin     ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
14276c6b9e74SPeter Brune     ierr = SNESSetFunction(ts->snes,PETSC_NULL,SNESTSFormFunction,ts);CHKERRQ(ierr);
1428d372ba47SLisandro Dalcin     ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr);
1429d372ba47SLisandro Dalcin     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
1430496e6a7aSJed Brown     if (ts->dm) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
14319e2a6581SJed Brown     if (ts->problem_type == TS_LINEAR) {
14329e2a6581SJed Brown       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
14339e2a6581SJed Brown     }
1434d372ba47SLisandro Dalcin   }
1435d763cef2SBarry Smith   *snes = ts->snes;
1436d763cef2SBarry Smith   PetscFunctionReturn(0);
1437d763cef2SBarry Smith }
1438d763cef2SBarry Smith 
14394a2ae208SSatish Balay #undef __FUNCT__
144094b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
1441d8e5e3e6SSatish Balay /*@
144294b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
1443d763cef2SBarry Smith    a TS (timestepper) context.
1444d763cef2SBarry Smith 
144594b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1446d763cef2SBarry Smith 
1447d763cef2SBarry Smith    Input Parameter:
1448d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1449d763cef2SBarry Smith 
1450d763cef2SBarry Smith    Output Parameter:
145194b7f48cSBarry Smith .  ksp - the nonlinear solver context
1452d763cef2SBarry Smith 
1453d763cef2SBarry Smith    Notes:
145494b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1455d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1456d763cef2SBarry Smith    KSP and PC contexts as well.
1457d763cef2SBarry Smith 
145894b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
145994b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1460d763cef2SBarry Smith 
1461d763cef2SBarry Smith    Level: beginner
1462d763cef2SBarry Smith 
146394b7f48cSBarry Smith .keywords: timestep, get, KSP
1464d763cef2SBarry Smith @*/
14657087cfbeSBarry Smith PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1466d763cef2SBarry Smith {
1467d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1468089b2837SJed Brown   SNES           snes;
1469d372ba47SLisandro Dalcin 
1470d763cef2SBarry Smith   PetscFunctionBegin;
14710700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
14724482741eSBarry Smith   PetscValidPointer(ksp,2);
147317186662SBarry Smith   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1474e32f2f54SBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1475089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1476089b2837SJed Brown   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
1477d763cef2SBarry Smith   PetscFunctionReturn(0);
1478d763cef2SBarry Smith }
1479d763cef2SBarry Smith 
1480d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1481d763cef2SBarry Smith 
14824a2ae208SSatish Balay #undef __FUNCT__
1483adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1484adb62b0dSMatthew Knepley /*@
1485adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1486adb62b0dSMatthew Knepley    maximum time for iteration.
1487adb62b0dSMatthew Knepley 
14883f9fe445SBarry Smith    Not Collective
1489adb62b0dSMatthew Knepley 
1490adb62b0dSMatthew Knepley    Input Parameters:
1491adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1492adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1493adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1494adb62b0dSMatthew Knepley 
1495adb62b0dSMatthew Knepley    Level: intermediate
1496adb62b0dSMatthew Knepley 
1497adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1498adb62b0dSMatthew Knepley @*/
14997087cfbeSBarry Smith PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1500adb62b0dSMatthew Knepley {
1501adb62b0dSMatthew Knepley   PetscFunctionBegin;
15020700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1503abc0a331SBarry Smith   if (maxsteps) {
15044482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1505adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1506adb62b0dSMatthew Knepley   }
1507abc0a331SBarry Smith   if (maxtime) {
15084482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1509adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1510adb62b0dSMatthew Knepley   }
1511adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1512adb62b0dSMatthew Knepley }
1513adb62b0dSMatthew Knepley 
1514adb62b0dSMatthew Knepley #undef __FUNCT__
15154a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1516d763cef2SBarry Smith /*@
1517d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1518d763cef2SBarry Smith    maximum time for iteration.
1519d763cef2SBarry Smith 
15203f9fe445SBarry Smith    Logically Collective on TS
1521d763cef2SBarry Smith 
1522d763cef2SBarry Smith    Input Parameters:
1523d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1524d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1525d763cef2SBarry Smith -  maxtime - final time to iterate to
1526d763cef2SBarry Smith 
1527d763cef2SBarry Smith    Options Database Keys:
1528d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
15293bca7d26SBarry Smith .  -ts_final_time <maxtime> - Sets maxtime
1530d763cef2SBarry Smith 
1531d763cef2SBarry Smith    Notes:
1532d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1533d763cef2SBarry Smith 
1534d763cef2SBarry Smith    Level: intermediate
1535d763cef2SBarry Smith 
1536d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1537a43b19c4SJed Brown 
1538a43b19c4SJed Brown .seealso: TSSetExactFinalTime()
1539d763cef2SBarry Smith @*/
15407087cfbeSBarry Smith PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1541d763cef2SBarry Smith {
1542d763cef2SBarry Smith   PetscFunctionBegin;
15430700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1544c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1545c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,maxtime,2);
154639b7ec4bSSean Farley   if (maxsteps >= 0) ts->max_steps = maxsteps;
154739b7ec4bSSean Farley   if (maxtime != PETSC_DEFAULT) ts->max_time  = maxtime;
1548d763cef2SBarry Smith   PetscFunctionReturn(0);
1549d763cef2SBarry Smith }
1550d763cef2SBarry Smith 
15514a2ae208SSatish Balay #undef __FUNCT__
15524a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1553d763cef2SBarry Smith /*@
1554d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1555d763cef2SBarry Smith    for use by the TS routines.
1556d763cef2SBarry Smith 
15573f9fe445SBarry Smith    Logically Collective on TS and Vec
1558d763cef2SBarry Smith 
1559d763cef2SBarry Smith    Input Parameters:
1560d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1561d763cef2SBarry Smith -  x - the solution vector
1562d763cef2SBarry Smith 
1563d763cef2SBarry Smith    Level: beginner
1564d763cef2SBarry Smith 
1565d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1566d763cef2SBarry Smith @*/
15677087cfbeSBarry Smith PetscErrorCode  TSSetSolution(TS ts,Vec x)
1568d763cef2SBarry Smith {
15698737fe31SLisandro Dalcin   PetscErrorCode ierr;
1570496e6a7aSJed Brown   DM             dm;
15718737fe31SLisandro Dalcin 
1572d763cef2SBarry Smith   PetscFunctionBegin;
15730700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
15740700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
15758737fe31SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
15766bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
15778737fe31SLisandro Dalcin   ts->vec_sol = x;
1578496e6a7aSJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1579496e6a7aSJed Brown   ierr = DMShellSetGlobalVector(dm,x);CHKERRQ(ierr);
1580d763cef2SBarry Smith   PetscFunctionReturn(0);
1581d763cef2SBarry Smith }
1582d763cef2SBarry Smith 
1583e74ef692SMatthew Knepley #undef __FUNCT__
1584e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1585ac226902SBarry Smith /*@C
1586000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
15873f2090d5SJed Brown   called once at the beginning of each time step.
1588000e7ae3SMatthew Knepley 
15893f9fe445SBarry Smith   Logically Collective on TS
1590000e7ae3SMatthew Knepley 
1591000e7ae3SMatthew Knepley   Input Parameters:
1592000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1593000e7ae3SMatthew Knepley - func - The function
1594000e7ae3SMatthew Knepley 
1595000e7ae3SMatthew Knepley   Calling sequence of func:
1596000e7ae3SMatthew Knepley . func (TS ts);
1597000e7ae3SMatthew Knepley 
1598000e7ae3SMatthew Knepley   Level: intermediate
1599000e7ae3SMatthew Knepley 
1600b8123daeSJed Brown   Note:
1601b8123daeSJed Brown   If a step is rejected, TSStep() will call this routine again before each attempt.
1602b8123daeSJed Brown   The last completed time step number can be queried using TSGetTimeStepNumber(), the
1603b8123daeSJed Brown   size of the step being attempted can be obtained using TSGetTimeStep().
1604b8123daeSJed Brown 
1605000e7ae3SMatthew Knepley .keywords: TS, timestep
1606b8123daeSJed Brown .seealso: TSSetPreStage(), TSSetPostStep(), TSStep()
1607000e7ae3SMatthew Knepley @*/
16087087cfbeSBarry Smith PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1609000e7ae3SMatthew Knepley {
1610000e7ae3SMatthew Knepley   PetscFunctionBegin;
16110700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1612000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1613000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1614000e7ae3SMatthew Knepley }
1615000e7ae3SMatthew Knepley 
1616e74ef692SMatthew Knepley #undef __FUNCT__
16173f2090d5SJed Brown #define __FUNCT__ "TSPreStep"
161809ee8438SJed Brown /*@
16193f2090d5SJed Brown   TSPreStep - Runs the user-defined pre-step function.
16203f2090d5SJed Brown 
16213f2090d5SJed Brown   Collective on TS
16223f2090d5SJed Brown 
16233f2090d5SJed Brown   Input Parameters:
16243f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
16253f2090d5SJed Brown 
16263f2090d5SJed Brown   Notes:
16273f2090d5SJed Brown   TSPreStep() is typically used within time stepping implementations,
16283f2090d5SJed Brown   so most users would not generally call this routine themselves.
16293f2090d5SJed Brown 
16303f2090d5SJed Brown   Level: developer
16313f2090d5SJed Brown 
16323f2090d5SJed Brown .keywords: TS, timestep
1633b8123daeSJed Brown .seealso: TSSetPreStep(), TSPreStage(), TSPostStep()
16343f2090d5SJed Brown @*/
16357087cfbeSBarry Smith PetscErrorCode  TSPreStep(TS ts)
16363f2090d5SJed Brown {
16373f2090d5SJed Brown   PetscErrorCode ierr;
16383f2090d5SJed Brown 
16393f2090d5SJed Brown   PetscFunctionBegin;
16400700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
164172ac3e02SJed Brown   if (ts->ops->prestep) {
16423f2090d5SJed Brown     PetscStackPush("TS PreStep function");
16433f2090d5SJed Brown     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
16443f2090d5SJed Brown     PetscStackPop;
1645312ce896SJed Brown   }
16463f2090d5SJed Brown   PetscFunctionReturn(0);
16473f2090d5SJed Brown }
16483f2090d5SJed Brown 
16493f2090d5SJed Brown #undef __FUNCT__
1650b8123daeSJed Brown #define __FUNCT__ "TSSetPreStage"
1651b8123daeSJed Brown /*@C
1652b8123daeSJed Brown   TSSetPreStage - Sets the general-purpose function
1653b8123daeSJed Brown   called once at the beginning of each stage.
1654b8123daeSJed Brown 
1655b8123daeSJed Brown   Logically Collective on TS
1656b8123daeSJed Brown 
1657b8123daeSJed Brown   Input Parameters:
1658b8123daeSJed Brown + ts   - The TS context obtained from TSCreate()
1659b8123daeSJed Brown - func - The function
1660b8123daeSJed Brown 
1661b8123daeSJed Brown   Calling sequence of func:
1662b8123daeSJed Brown . PetscErrorCode func(TS ts, PetscReal stagetime);
1663b8123daeSJed Brown 
1664b8123daeSJed Brown   Level: intermediate
1665b8123daeSJed Brown 
1666b8123daeSJed Brown   Note:
1667b8123daeSJed Brown   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
1668b8123daeSJed Brown   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
1669b8123daeSJed Brown   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
1670b8123daeSJed Brown 
1671b8123daeSJed Brown .keywords: TS, timestep
1672b8123daeSJed Brown .seealso: TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
1673b8123daeSJed Brown @*/
1674b8123daeSJed Brown PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
1675b8123daeSJed Brown {
1676b8123daeSJed Brown   PetscFunctionBegin;
1677b8123daeSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1678b8123daeSJed Brown   ts->ops->prestage = func;
1679b8123daeSJed Brown   PetscFunctionReturn(0);
1680b8123daeSJed Brown }
1681b8123daeSJed Brown 
1682b8123daeSJed Brown #undef __FUNCT__
1683b8123daeSJed Brown #define __FUNCT__ "TSPreStage"
1684b8123daeSJed Brown /*@
1685b8123daeSJed Brown   TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()
1686b8123daeSJed Brown 
1687b8123daeSJed Brown   Collective on TS
1688b8123daeSJed Brown 
1689b8123daeSJed Brown   Input Parameters:
1690b8123daeSJed Brown . ts   - The TS context obtained from TSCreate()
1691b8123daeSJed Brown 
1692b8123daeSJed Brown   Notes:
1693b8123daeSJed Brown   TSPreStage() is typically used within time stepping implementations,
1694b8123daeSJed Brown   most users would not generally call this routine themselves.
1695b8123daeSJed Brown 
1696b8123daeSJed Brown   Level: developer
1697b8123daeSJed Brown 
1698b8123daeSJed Brown .keywords: TS, timestep
1699b8123daeSJed Brown .seealso: TSSetPreStep(), TSPreStep(), TSPostStep()
1700b8123daeSJed Brown @*/
1701b8123daeSJed Brown PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
1702b8123daeSJed Brown {
1703b8123daeSJed Brown   PetscErrorCode ierr;
1704b8123daeSJed Brown 
1705b8123daeSJed Brown   PetscFunctionBegin;
1706b8123daeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1707b8123daeSJed Brown   if (ts->ops->prestage) {
1708b8123daeSJed Brown     PetscStackPush("TS PreStage function");
1709b8123daeSJed Brown     ierr = (*ts->ops->prestage)(ts,stagetime);CHKERRQ(ierr);
1710b8123daeSJed Brown     PetscStackPop;
1711b8123daeSJed Brown   }
1712b8123daeSJed Brown   PetscFunctionReturn(0);
1713b8123daeSJed Brown }
1714b8123daeSJed Brown 
1715b8123daeSJed Brown #undef __FUNCT__
1716e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1717ac226902SBarry Smith /*@C
1718000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
17193f2090d5SJed Brown   called once at the end of each time step.
1720000e7ae3SMatthew Knepley 
17213f9fe445SBarry Smith   Logically Collective on TS
1722000e7ae3SMatthew Knepley 
1723000e7ae3SMatthew Knepley   Input Parameters:
1724000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1725000e7ae3SMatthew Knepley - func - The function
1726000e7ae3SMatthew Knepley 
1727000e7ae3SMatthew Knepley   Calling sequence of func:
1728b8123daeSJed Brown $ func (TS ts);
1729000e7ae3SMatthew Knepley 
1730000e7ae3SMatthew Knepley   Level: intermediate
1731000e7ae3SMatthew Knepley 
1732000e7ae3SMatthew Knepley .keywords: TS, timestep
1733b8123daeSJed Brown .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime()
1734000e7ae3SMatthew Knepley @*/
17357087cfbeSBarry Smith PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1736000e7ae3SMatthew Knepley {
1737000e7ae3SMatthew Knepley   PetscFunctionBegin;
17380700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1739000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1740000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1741000e7ae3SMatthew Knepley }
1742000e7ae3SMatthew Knepley 
1743e74ef692SMatthew Knepley #undef __FUNCT__
17443f2090d5SJed Brown #define __FUNCT__ "TSPostStep"
174509ee8438SJed Brown /*@
17463f2090d5SJed Brown   TSPostStep - Runs the user-defined post-step function.
17473f2090d5SJed Brown 
17483f2090d5SJed Brown   Collective on TS
17493f2090d5SJed Brown 
17503f2090d5SJed Brown   Input Parameters:
17513f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
17523f2090d5SJed Brown 
17533f2090d5SJed Brown   Notes:
17543f2090d5SJed Brown   TSPostStep() is typically used within time stepping implementations,
17553f2090d5SJed Brown   so most users would not generally call this routine themselves.
17563f2090d5SJed Brown 
17573f2090d5SJed Brown   Level: developer
17583f2090d5SJed Brown 
17593f2090d5SJed Brown .keywords: TS, timestep
17603f2090d5SJed Brown @*/
17617087cfbeSBarry Smith PetscErrorCode  TSPostStep(TS ts)
17623f2090d5SJed Brown {
17633f2090d5SJed Brown   PetscErrorCode ierr;
17643f2090d5SJed Brown 
17653f2090d5SJed Brown   PetscFunctionBegin;
17660700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
176772ac3e02SJed Brown   if (ts->ops->poststep) {
17683f2090d5SJed Brown     PetscStackPush("TS PostStep function");
17693f2090d5SJed Brown     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
17703f2090d5SJed Brown     PetscStackPop;
177172ac3e02SJed Brown   }
17723f2090d5SJed Brown   PetscFunctionReturn(0);
17733f2090d5SJed Brown }
17743f2090d5SJed Brown 
1775d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1776d763cef2SBarry Smith 
17774a2ae208SSatish Balay #undef __FUNCT__
1778a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet"
1779d763cef2SBarry Smith /*@C
1780a6570f20SBarry Smith    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1781d763cef2SBarry Smith    timestep to display the iteration's  progress.
1782d763cef2SBarry Smith 
17833f9fe445SBarry Smith    Logically Collective on TS
1784d763cef2SBarry Smith 
1785d763cef2SBarry Smith    Input Parameters:
1786d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1787e213d8f1SJed Brown .  monitor - monitoring routine
1788329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1789b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1790b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1791b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1792d763cef2SBarry Smith 
1793e213d8f1SJed Brown    Calling sequence of monitor:
1794e213d8f1SJed Brown $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1795d763cef2SBarry Smith 
1796d763cef2SBarry Smith +    ts - the TS context
1797d763cef2SBarry Smith .    steps - iteration number
17981f06c33eSBarry Smith .    time - current time
1799d763cef2SBarry Smith .    x - current iterate
1800d763cef2SBarry Smith -    mctx - [optional] monitoring context
1801d763cef2SBarry Smith 
1802d763cef2SBarry Smith    Notes:
1803d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1804d763cef2SBarry Smith    already has been loaded.
1805d763cef2SBarry Smith 
1806025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each TS object
1807025f1a04SBarry Smith 
1808d763cef2SBarry Smith    Level: intermediate
1809d763cef2SBarry Smith 
1810d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1811d763cef2SBarry Smith 
1812a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel()
1813d763cef2SBarry Smith @*/
1814c2efdce3SBarry Smith PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1815d763cef2SBarry Smith {
1816d763cef2SBarry Smith   PetscFunctionBegin;
18170700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
181817186662SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1819d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1820329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1821d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1822d763cef2SBarry Smith   PetscFunctionReturn(0);
1823d763cef2SBarry Smith }
1824d763cef2SBarry Smith 
18254a2ae208SSatish Balay #undef __FUNCT__
1826a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel"
1827d763cef2SBarry Smith /*@C
1828a6570f20SBarry Smith    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1829d763cef2SBarry Smith 
18303f9fe445SBarry Smith    Logically Collective on TS
1831d763cef2SBarry Smith 
1832d763cef2SBarry Smith    Input Parameters:
1833d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1834d763cef2SBarry Smith 
1835d763cef2SBarry Smith    Notes:
1836d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1837d763cef2SBarry Smith 
1838d763cef2SBarry Smith    Level: intermediate
1839d763cef2SBarry Smith 
1840d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1841d763cef2SBarry Smith 
1842a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet()
1843d763cef2SBarry Smith @*/
18447087cfbeSBarry Smith PetscErrorCode  TSMonitorCancel(TS ts)
1845d763cef2SBarry Smith {
1846d952e501SBarry Smith   PetscErrorCode ierr;
1847d952e501SBarry Smith   PetscInt       i;
1848d952e501SBarry Smith 
1849d763cef2SBarry Smith   PetscFunctionBegin;
18500700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1851d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1852d952e501SBarry Smith     if (ts->mdestroy[i]) {
18533c4aec1bSBarry Smith       ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
1854d952e501SBarry Smith     }
1855d952e501SBarry Smith   }
1856d763cef2SBarry Smith   ts->numbermonitors = 0;
1857d763cef2SBarry Smith   PetscFunctionReturn(0);
1858d763cef2SBarry Smith }
1859d763cef2SBarry Smith 
18604a2ae208SSatish Balay #undef __FUNCT__
1861a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault"
1862d8e5e3e6SSatish Balay /*@
1863a6570f20SBarry Smith    TSMonitorDefault - Sets the Default monitor
18645516499fSSatish Balay 
18655516499fSSatish Balay    Level: intermediate
186641251cbbSSatish Balay 
18675516499fSSatish Balay .keywords: TS, set, monitor
18685516499fSSatish Balay 
186941251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet()
187041251cbbSSatish Balay @*/
1871649052a6SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1872d763cef2SBarry Smith {
1873dfbe8321SBarry Smith   PetscErrorCode ierr;
1874649052a6SBarry Smith   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1875d132466eSBarry Smith 
1876d763cef2SBarry Smith   PetscFunctionBegin;
1877649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1878971832b6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr);
1879649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1880d763cef2SBarry Smith   PetscFunctionReturn(0);
1881d763cef2SBarry Smith }
1882d763cef2SBarry Smith 
18834a2ae208SSatish Balay #undef __FUNCT__
1884cd652676SJed Brown #define __FUNCT__ "TSSetRetainStages"
1885cd652676SJed Brown /*@
1886cd652676SJed Brown    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
1887cd652676SJed Brown 
1888cd652676SJed Brown    Logically Collective on TS
1889cd652676SJed Brown 
1890cd652676SJed Brown    Input Argument:
1891cd652676SJed Brown .  ts - time stepping context
1892cd652676SJed Brown 
1893cd652676SJed Brown    Output Argument:
1894cd652676SJed Brown .  flg - PETSC_TRUE or PETSC_FALSE
1895cd652676SJed Brown 
1896cd652676SJed Brown    Level: intermediate
1897cd652676SJed Brown 
1898cd652676SJed Brown .keywords: TS, set
1899cd652676SJed Brown 
1900cd652676SJed Brown .seealso: TSInterpolate(), TSSetPostStep()
1901cd652676SJed Brown @*/
1902cd652676SJed Brown PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
1903cd652676SJed Brown {
1904cd652676SJed Brown 
1905cd652676SJed Brown   PetscFunctionBegin;
1906cd652676SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1907cd652676SJed Brown   ts->retain_stages = flg;
1908cd652676SJed Brown   PetscFunctionReturn(0);
1909cd652676SJed Brown }
1910cd652676SJed Brown 
1911cd652676SJed Brown #undef __FUNCT__
1912cd652676SJed Brown #define __FUNCT__ "TSInterpolate"
1913cd652676SJed Brown /*@
1914cd652676SJed Brown    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
1915cd652676SJed Brown 
1916cd652676SJed Brown    Collective on TS
1917cd652676SJed Brown 
1918cd652676SJed Brown    Input Argument:
1919cd652676SJed Brown +  ts - time stepping context
1920cd652676SJed Brown -  t - time to interpolate to
1921cd652676SJed Brown 
1922cd652676SJed Brown    Output Argument:
1923cd652676SJed Brown .  X - state at given time
1924cd652676SJed Brown 
1925cd652676SJed Brown    Notes:
1926cd652676SJed Brown    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
1927cd652676SJed Brown 
1928cd652676SJed Brown    Level: intermediate
1929cd652676SJed Brown 
1930cd652676SJed Brown    Developer Notes:
1931cd652676SJed Brown    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
1932cd652676SJed Brown 
1933cd652676SJed Brown .keywords: TS, set
1934cd652676SJed Brown 
1935cd652676SJed Brown .seealso: TSSetRetainStages(), TSSetPostStep()
1936cd652676SJed Brown @*/
1937cd652676SJed Brown PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
1938cd652676SJed Brown {
1939cd652676SJed Brown   PetscErrorCode ierr;
1940cd652676SJed Brown 
1941cd652676SJed Brown   PetscFunctionBegin;
1942cd652676SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1943d8cd7023SBarry 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);
1944cd652676SJed Brown   if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
1945cd652676SJed Brown   ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr);
1946cd652676SJed Brown   PetscFunctionReturn(0);
1947cd652676SJed Brown }
1948cd652676SJed Brown 
1949cd652676SJed Brown #undef __FUNCT__
19504a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1951d763cef2SBarry Smith /*@
19526d9e5789SSean Farley    TSStep - Steps one time step
1953d763cef2SBarry Smith 
1954d763cef2SBarry Smith    Collective on TS
1955d763cef2SBarry Smith 
1956d763cef2SBarry Smith    Input Parameter:
1957d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1958d763cef2SBarry Smith 
19596d9e5789SSean Farley    Level: intermediate
1960d763cef2SBarry Smith 
1961b8123daeSJed Brown    Notes:
1962b8123daeSJed Brown    The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
1963b8123daeSJed Brown    be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.
1964b8123daeSJed Brown 
1965d763cef2SBarry Smith .keywords: TS, timestep, solve
1966d763cef2SBarry Smith 
1967b8123daeSJed Brown .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage()
1968d763cef2SBarry Smith @*/
1969193ac0bcSJed Brown PetscErrorCode  TSStep(TS ts)
1970d763cef2SBarry Smith {
1971362cd11cSLisandro Dalcin   PetscReal      ptime_prev;
1972dfbe8321SBarry Smith   PetscErrorCode ierr;
1973d763cef2SBarry Smith 
1974d763cef2SBarry Smith   PetscFunctionBegin;
19750700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1976d405a339SMatthew Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
1977d405a339SMatthew Knepley 
1978362cd11cSLisandro Dalcin   ts->reason = TS_CONVERGED_ITERATING;
1979362cd11cSLisandro Dalcin 
1980362cd11cSLisandro Dalcin   ptime_prev = ts->ptime;
1981d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1982193ac0bcSJed Brown   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
1983d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1984362cd11cSLisandro Dalcin   ts->time_step_prev = ts->ptime - ptime_prev;
1985362cd11cSLisandro Dalcin 
1986362cd11cSLisandro Dalcin   if (ts->reason < 0) {
1987cef5090cSJed Brown     if (ts->errorifstepfailed) {
1988cef5090cSJed Brown       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
1989cef5090cSJed Brown         SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
1990cef5090cSJed Brown       } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
1991cef5090cSJed Brown     }
1992362cd11cSLisandro Dalcin   } else if (!ts->reason) {
1993362cd11cSLisandro Dalcin     if (ts->steps >= ts->max_steps)
1994362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_ITS;
1995362cd11cSLisandro Dalcin     else if (ts->ptime >= ts->max_time)
1996362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_TIME;
1997362cd11cSLisandro Dalcin   }
1998362cd11cSLisandro Dalcin 
1999d763cef2SBarry Smith   PetscFunctionReturn(0);
2000d763cef2SBarry Smith }
2001d763cef2SBarry Smith 
20024a2ae208SSatish Balay #undef __FUNCT__
200305175c85SJed Brown #define __FUNCT__ "TSEvaluateStep"
200405175c85SJed Brown /*@
200505175c85SJed Brown    TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
200605175c85SJed Brown 
20071c3436cfSJed Brown    Collective on TS
200805175c85SJed Brown 
200905175c85SJed Brown    Input Arguments:
20101c3436cfSJed Brown +  ts - time stepping context
20111c3436cfSJed Brown .  order - desired order of accuracy
20121c3436cfSJed Brown -  done - whether the step was evaluated at this order (pass PETSC_NULL to generate an error if not available)
201305175c85SJed Brown 
201405175c85SJed Brown    Output Arguments:
20151c3436cfSJed Brown .  X - state at the end of the current step
201605175c85SJed Brown 
201705175c85SJed Brown    Level: advanced
201805175c85SJed Brown 
2019108c343cSJed Brown    Notes:
2020108c343cSJed Brown    This function cannot be called until all stages have been evaluated.
2021108c343cSJed Brown    It is normally called by adaptive controllers before a step has been accepted and may also be called by the user after TSStep() has returned.
2022108c343cSJed Brown 
20231c3436cfSJed Brown .seealso: TSStep(), TSAdapt
202405175c85SJed Brown @*/
20251c3436cfSJed Brown PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done)
202605175c85SJed Brown {
202705175c85SJed Brown   PetscErrorCode ierr;
202805175c85SJed Brown 
202905175c85SJed Brown   PetscFunctionBegin;
203005175c85SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
203105175c85SJed Brown   PetscValidType(ts,1);
203205175c85SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
20331c3436cfSJed Brown   if (!ts->ops->evaluatestep) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
20341c3436cfSJed Brown   ierr = (*ts->ops->evaluatestep)(ts,order,X,done);CHKERRQ(ierr);
203505175c85SJed Brown   PetscFunctionReturn(0);
203605175c85SJed Brown }
203705175c85SJed Brown 
203805175c85SJed Brown #undef __FUNCT__
20396a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve"
20406a4d4014SLisandro Dalcin /*@
20416a4d4014SLisandro Dalcin    TSSolve - Steps the requested number of timesteps.
20426a4d4014SLisandro Dalcin 
20436a4d4014SLisandro Dalcin    Collective on TS
20446a4d4014SLisandro Dalcin 
20456a4d4014SLisandro Dalcin    Input Parameter:
20466a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
2047362cd11cSLisandro Dalcin -  x - the solution vector
20486a4d4014SLisandro Dalcin 
20495a3a76d0SJed Brown    Output Parameter:
20505a3a76d0SJed Brown .  ftime - time of the state vector x upon completion
20515a3a76d0SJed Brown 
20526a4d4014SLisandro Dalcin    Level: beginner
20536a4d4014SLisandro Dalcin 
20545a3a76d0SJed Brown    Notes:
20555a3a76d0SJed Brown    The final time returned by this function may be different from the time of the internally
20565a3a76d0SJed Brown    held state accessible by TSGetSolution() and TSGetTime() because the method may have
20575a3a76d0SJed Brown    stepped over the final time.
20585a3a76d0SJed Brown 
20596a4d4014SLisandro Dalcin .keywords: TS, timestep, solve
20606a4d4014SLisandro Dalcin 
20616a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep()
20626a4d4014SLisandro Dalcin @*/
20635a3a76d0SJed Brown PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime)
20646a4d4014SLisandro Dalcin {
20654d7d938eSLisandro Dalcin   PetscBool      flg;
20664d7d938eSLisandro Dalcin   char           filename[PETSC_MAX_PATH_LEN];
20674d7d938eSLisandro Dalcin   PetscViewer    viewer;
20686a4d4014SLisandro Dalcin   PetscErrorCode ierr;
2069f22f69f0SBarry Smith 
20706a4d4014SLisandro Dalcin   PetscFunctionBegin;
20710700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2072193ac0bcSJed Brown   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
20735a3a76d0SJed Brown   if (ts->exact_final_time) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
20745a3a76d0SJed Brown     if (!ts->vec_sol || x == ts->vec_sol) {
20755a3a76d0SJed Brown       Vec y;
20765a3a76d0SJed Brown       ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
20775a3a76d0SJed Brown       ierr = TSSetSolution(ts,y);CHKERRQ(ierr);
20785a3a76d0SJed Brown       ierr = VecDestroy(&y);CHKERRQ(ierr); /* grant ownership */
20795a3a76d0SJed Brown     }
2080362cd11cSLisandro Dalcin     ierr = VecCopy(x,ts->vec_sol);CHKERRQ(ierr);
20815a3a76d0SJed Brown   } else {
2082193ac0bcSJed Brown     ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
20835a3a76d0SJed Brown   }
2084b5d403baSSean Farley   ierr = TSSetUp(ts);CHKERRQ(ierr);
20856a4d4014SLisandro Dalcin   /* reset time step and iteration counters */
2086193ac0bcSJed Brown   ts->steps = 0;
20875ef26d82SJed Brown   ts->ksp_its = 0;
20885ef26d82SJed Brown   ts->snes_its = 0;
2089c610991cSLisandro Dalcin   ts->num_snes_failures = 0;
2090c610991cSLisandro Dalcin   ts->reject = 0;
2091193ac0bcSJed Brown   ts->reason = TS_CONVERGED_ITERATING;
2092193ac0bcSJed Brown 
2093193ac0bcSJed Brown   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
2094193ac0bcSJed Brown     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
20955a3a76d0SJed Brown     ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
2096a5b82c69SSean Farley     if (ftime) *ftime = ts->ptime;
2097193ac0bcSJed Brown   } else {
20986a4d4014SLisandro Dalcin     /* steps the requested number of timesteps. */
2099362cd11cSLisandro Dalcin     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2100362cd11cSLisandro Dalcin     if (ts->steps >= ts->max_steps)
2101362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_ITS;
2102362cd11cSLisandro Dalcin     else if (ts->ptime >= ts->max_time)
2103362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_TIME;
2104e1a7a14fSJed Brown     while (!ts->reason) {
2105193ac0bcSJed Brown       ierr = TSStep(ts);CHKERRQ(ierr);
2106193ac0bcSJed Brown       ierr = TSPostStep(ts);CHKERRQ(ierr);
2107193ac0bcSJed Brown       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2108193ac0bcSJed Brown     }
2109362cd11cSLisandro Dalcin     if (ts->exact_final_time && ts->ptime > ts->max_time) {
2110a43b19c4SJed Brown       ierr = TSInterpolate(ts,ts->max_time,x);CHKERRQ(ierr);
21115a3a76d0SJed Brown       if (ftime) *ftime = ts->max_time;
21120574a7fbSJed Brown     } else {
21130574a7fbSJed Brown       ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
21140574a7fbSJed Brown       if (ftime) *ftime = ts->ptime;
21150574a7fbSJed Brown     }
2116193ac0bcSJed Brown   }
21174d7d938eSLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
21184d7d938eSLisandro Dalcin   if (flg && !PetscPreLoadingOn) {
21194d7d938eSLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr);
21204d7d938eSLisandro Dalcin     ierr = TSView(ts,viewer);CHKERRQ(ierr);
21214d7d938eSLisandro Dalcin     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2122193ac0bcSJed Brown   }
21236a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
21246a4d4014SLisandro Dalcin }
21256a4d4014SLisandro Dalcin 
21266a4d4014SLisandro Dalcin #undef __FUNCT__
21274a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
2128228d79bcSJed Brown /*@
2129228d79bcSJed Brown    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
2130228d79bcSJed Brown 
2131228d79bcSJed Brown    Collective on TS
2132228d79bcSJed Brown 
2133228d79bcSJed Brown    Input Parameters:
2134228d79bcSJed Brown +  ts - time stepping context obtained from TSCreate()
2135228d79bcSJed Brown .  step - step number that has just completed
2136228d79bcSJed Brown .  ptime - model time of the state
2137228d79bcSJed Brown -  x - state at the current model time
2138228d79bcSJed Brown 
2139228d79bcSJed Brown    Notes:
2140228d79bcSJed Brown    TSMonitor() is typically used within the time stepping implementations.
2141228d79bcSJed Brown    Users might call this function when using the TSStep() interface instead of TSSolve().
2142228d79bcSJed Brown 
2143228d79bcSJed Brown    Level: advanced
2144228d79bcSJed Brown 
2145228d79bcSJed Brown .keywords: TS, timestep
2146228d79bcSJed Brown @*/
2147a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
2148d763cef2SBarry Smith {
21496849ba73SBarry Smith   PetscErrorCode ierr;
2150a7cc72afSBarry Smith   PetscInt       i,n = ts->numbermonitors;
2151d763cef2SBarry Smith 
2152d763cef2SBarry Smith   PetscFunctionBegin;
2153d763cef2SBarry Smith   for (i=0; i<n; i++) {
2154142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
2155d763cef2SBarry Smith   }
2156d763cef2SBarry Smith   PetscFunctionReturn(0);
2157d763cef2SBarry Smith }
2158d763cef2SBarry Smith 
2159d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
2160d763cef2SBarry Smith 
21614a2ae208SSatish Balay #undef __FUNCT__
2162a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate"
2163d763cef2SBarry Smith /*@C
2164a6570f20SBarry Smith    TSMonitorLGCreate - Creates a line graph context for use with
2165d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
2166d763cef2SBarry Smith 
2167d763cef2SBarry Smith    Collective on TS
2168d763cef2SBarry Smith 
2169d763cef2SBarry Smith    Input Parameters:
2170d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
2171d763cef2SBarry Smith .  label - the title to put in the title bar
21727c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
2173d763cef2SBarry Smith -  m, n - the screen width and height in pixels
2174d763cef2SBarry Smith 
2175d763cef2SBarry Smith    Output Parameter:
2176d763cef2SBarry Smith .  draw - the drawing context
2177d763cef2SBarry Smith 
2178d763cef2SBarry Smith    Options Database Key:
2179a6570f20SBarry Smith .  -ts_monitor_draw - automatically sets line graph monitor
2180d763cef2SBarry Smith 
2181d763cef2SBarry Smith    Notes:
2182a6570f20SBarry Smith    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
2183d763cef2SBarry Smith 
2184d763cef2SBarry Smith    Level: intermediate
2185d763cef2SBarry Smith 
21867c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
2187d763cef2SBarry Smith 
2188a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet()
21897c922b88SBarry Smith 
2190d763cef2SBarry Smith @*/
21917087cfbeSBarry Smith PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2192d763cef2SBarry Smith {
2193b0a32e0cSBarry Smith   PetscDraw      win;
2194dfbe8321SBarry Smith   PetscErrorCode ierr;
2195d763cef2SBarry Smith 
2196d763cef2SBarry Smith   PetscFunctionBegin;
2197b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
2198b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
2199b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
2200b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
2201d763cef2SBarry Smith 
220252e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
2203d763cef2SBarry Smith   PetscFunctionReturn(0);
2204d763cef2SBarry Smith }
2205d763cef2SBarry Smith 
22064a2ae208SSatish Balay #undef __FUNCT__
2207a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG"
2208a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
2209d763cef2SBarry Smith {
2210b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
221187828ca2SBarry Smith   PetscReal      x,y = ptime;
2212dfbe8321SBarry Smith   PetscErrorCode ierr;
2213d763cef2SBarry Smith 
2214d763cef2SBarry Smith   PetscFunctionBegin;
22157c922b88SBarry Smith   if (!monctx) {
22167c922b88SBarry Smith     MPI_Comm    comm;
2217b0a32e0cSBarry Smith     PetscViewer viewer;
22187c922b88SBarry Smith 
22197c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
2220b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
2221b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
22227c922b88SBarry Smith   }
22237c922b88SBarry Smith 
2224b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
222587828ca2SBarry Smith   x = (PetscReal)n;
2226b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2227d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
2228b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2229d763cef2SBarry Smith   }
2230d763cef2SBarry Smith   PetscFunctionReturn(0);
2231d763cef2SBarry Smith }
2232d763cef2SBarry Smith 
22334a2ae208SSatish Balay #undef __FUNCT__
2234a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy"
2235d763cef2SBarry Smith /*@C
2236a6570f20SBarry Smith    TSMonitorLGDestroy - Destroys a line graph context that was created
2237a6570f20SBarry Smith    with TSMonitorLGCreate().
2238d763cef2SBarry Smith 
2239b0a32e0cSBarry Smith    Collective on PetscDrawLG
2240d763cef2SBarry Smith 
2241d763cef2SBarry Smith    Input Parameter:
2242d763cef2SBarry Smith .  draw - the drawing context
2243d763cef2SBarry Smith 
2244d763cef2SBarry Smith    Level: intermediate
2245d763cef2SBarry Smith 
2246d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
2247d763cef2SBarry Smith 
2248a6570f20SBarry Smith .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
2249d763cef2SBarry Smith @*/
22506bf464f9SBarry Smith PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
2251d763cef2SBarry Smith {
2252b0a32e0cSBarry Smith   PetscDraw      draw;
2253dfbe8321SBarry Smith   PetscErrorCode ierr;
2254d763cef2SBarry Smith 
2255d763cef2SBarry Smith   PetscFunctionBegin;
22566bf464f9SBarry Smith   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
22576bf464f9SBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
2258b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
2259d763cef2SBarry Smith   PetscFunctionReturn(0);
2260d763cef2SBarry Smith }
2261d763cef2SBarry Smith 
22624a2ae208SSatish Balay #undef __FUNCT__
22634a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
2264d763cef2SBarry Smith /*@
2265b8123daeSJed Brown    TSGetTime - Gets the time of the most recently completed step.
2266d763cef2SBarry Smith 
2267d763cef2SBarry Smith    Not Collective
2268d763cef2SBarry Smith 
2269d763cef2SBarry Smith    Input Parameter:
2270d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
2271d763cef2SBarry Smith 
2272d763cef2SBarry Smith    Output Parameter:
2273d763cef2SBarry Smith .  t  - the current time
2274d763cef2SBarry Smith 
2275d763cef2SBarry Smith    Level: beginner
2276d763cef2SBarry Smith 
2277b8123daeSJed Brown    Note:
2278b8123daeSJed Brown    When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
2279b8123daeSJed Brown    TSSetPreStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.
2280b8123daeSJed Brown 
2281d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
2282d763cef2SBarry Smith 
2283d763cef2SBarry Smith .keywords: TS, get, time
2284d763cef2SBarry Smith @*/
22857087cfbeSBarry Smith PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
2286d763cef2SBarry Smith {
2287d763cef2SBarry Smith   PetscFunctionBegin;
22880700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2289f7cf8827SBarry Smith   PetscValidRealPointer(t,2);
2290d763cef2SBarry Smith   *t = ts->ptime;
2291d763cef2SBarry Smith   PetscFunctionReturn(0);
2292d763cef2SBarry Smith }
2293d763cef2SBarry Smith 
22944a2ae208SSatish Balay #undef __FUNCT__
22956a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime"
22966a4d4014SLisandro Dalcin /*@
22976a4d4014SLisandro Dalcin    TSSetTime - Allows one to reset the time.
22986a4d4014SLisandro Dalcin 
22993f9fe445SBarry Smith    Logically Collective on TS
23006a4d4014SLisandro Dalcin 
23016a4d4014SLisandro Dalcin    Input Parameters:
23026a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
23036a4d4014SLisandro Dalcin -  time - the time
23046a4d4014SLisandro Dalcin 
23056a4d4014SLisandro Dalcin    Level: intermediate
23066a4d4014SLisandro Dalcin 
23076a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration()
23086a4d4014SLisandro Dalcin 
23096a4d4014SLisandro Dalcin .keywords: TS, set, time
23106a4d4014SLisandro Dalcin @*/
23117087cfbeSBarry Smith PetscErrorCode  TSSetTime(TS ts, PetscReal t)
23126a4d4014SLisandro Dalcin {
23136a4d4014SLisandro Dalcin   PetscFunctionBegin;
23140700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2315c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,t,2);
23166a4d4014SLisandro Dalcin   ts->ptime = t;
23176a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
23186a4d4014SLisandro Dalcin }
23196a4d4014SLisandro Dalcin 
23206a4d4014SLisandro Dalcin #undef __FUNCT__
23214a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
2322d763cef2SBarry Smith /*@C
2323d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
2324d763cef2SBarry Smith    TS options in the database.
2325d763cef2SBarry Smith 
23263f9fe445SBarry Smith    Logically Collective on TS
2327d763cef2SBarry Smith 
2328d763cef2SBarry Smith    Input Parameter:
2329d763cef2SBarry Smith +  ts     - The TS context
2330d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
2331d763cef2SBarry Smith 
2332d763cef2SBarry Smith    Notes:
2333d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
2334d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
2335d763cef2SBarry Smith    hyphen.
2336d763cef2SBarry Smith 
2337d763cef2SBarry Smith    Level: advanced
2338d763cef2SBarry Smith 
2339d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
2340d763cef2SBarry Smith 
2341d763cef2SBarry Smith .seealso: TSSetFromOptions()
2342d763cef2SBarry Smith 
2343d763cef2SBarry Smith @*/
23447087cfbeSBarry Smith PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2345d763cef2SBarry Smith {
2346dfbe8321SBarry Smith   PetscErrorCode ierr;
2347089b2837SJed Brown   SNES           snes;
2348d763cef2SBarry Smith 
2349d763cef2SBarry Smith   PetscFunctionBegin;
23500700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2351d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2352089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2353089b2837SJed Brown   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2354d763cef2SBarry Smith   PetscFunctionReturn(0);
2355d763cef2SBarry Smith }
2356d763cef2SBarry Smith 
2357d763cef2SBarry Smith 
23584a2ae208SSatish Balay #undef __FUNCT__
23594a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
2360d763cef2SBarry Smith /*@C
2361d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2362d763cef2SBarry Smith    TS options in the database.
2363d763cef2SBarry Smith 
23643f9fe445SBarry Smith    Logically Collective on TS
2365d763cef2SBarry Smith 
2366d763cef2SBarry Smith    Input Parameter:
2367d763cef2SBarry Smith +  ts     - The TS context
2368d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
2369d763cef2SBarry Smith 
2370d763cef2SBarry Smith    Notes:
2371d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
2372d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
2373d763cef2SBarry Smith    hyphen.
2374d763cef2SBarry Smith 
2375d763cef2SBarry Smith    Level: advanced
2376d763cef2SBarry Smith 
2377d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
2378d763cef2SBarry Smith 
2379d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
2380d763cef2SBarry Smith 
2381d763cef2SBarry Smith @*/
23827087cfbeSBarry Smith PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2383d763cef2SBarry Smith {
2384dfbe8321SBarry Smith   PetscErrorCode ierr;
2385089b2837SJed Brown   SNES           snes;
2386d763cef2SBarry Smith 
2387d763cef2SBarry Smith   PetscFunctionBegin;
23880700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2389d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2390089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2391089b2837SJed Brown   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2392d763cef2SBarry Smith   PetscFunctionReturn(0);
2393d763cef2SBarry Smith }
2394d763cef2SBarry Smith 
23954a2ae208SSatish Balay #undef __FUNCT__
23964a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
2397d763cef2SBarry Smith /*@C
2398d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
2399d763cef2SBarry Smith    TS options in the database.
2400d763cef2SBarry Smith 
2401d763cef2SBarry Smith    Not Collective
2402d763cef2SBarry Smith 
2403d763cef2SBarry Smith    Input Parameter:
2404d763cef2SBarry Smith .  ts - The TS context
2405d763cef2SBarry Smith 
2406d763cef2SBarry Smith    Output Parameter:
2407d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
2408d763cef2SBarry Smith 
2409d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
2410d763cef2SBarry Smith    sufficient length to hold the prefix.
2411d763cef2SBarry Smith 
2412d763cef2SBarry Smith    Level: intermediate
2413d763cef2SBarry Smith 
2414d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
2415d763cef2SBarry Smith 
2416d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
2417d763cef2SBarry Smith @*/
24187087cfbeSBarry Smith PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2419d763cef2SBarry Smith {
2420dfbe8321SBarry Smith   PetscErrorCode ierr;
2421d763cef2SBarry Smith 
2422d763cef2SBarry Smith   PetscFunctionBegin;
24230700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
24244482741eSBarry Smith   PetscValidPointer(prefix,2);
2425d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2426d763cef2SBarry Smith   PetscFunctionReturn(0);
2427d763cef2SBarry Smith }
2428d763cef2SBarry Smith 
24294a2ae208SSatish Balay #undef __FUNCT__
24304a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
2431d763cef2SBarry Smith /*@C
2432d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2433d763cef2SBarry Smith 
2434d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
2435d763cef2SBarry Smith 
2436d763cef2SBarry Smith    Input Parameter:
2437d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
2438d763cef2SBarry Smith 
2439d763cef2SBarry Smith    Output Parameters:
2440d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
2441d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
2442089b2837SJed Brown .  func - Function to compute the Jacobian of the RHS
2443d763cef2SBarry Smith -  ctx - User-defined context for Jacobian evaluation routine
2444d763cef2SBarry Smith 
2445d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2446d763cef2SBarry Smith 
2447d763cef2SBarry Smith    Level: intermediate
2448d763cef2SBarry Smith 
244926d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2450d763cef2SBarry Smith 
2451d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
2452d763cef2SBarry Smith @*/
2453089b2837SJed Brown PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2454d763cef2SBarry Smith {
2455089b2837SJed Brown   PetscErrorCode ierr;
2456089b2837SJed Brown   SNES           snes;
245724989b8cSPeter Brune   DM             dm;
2458089b2837SJed Brown 
2459d763cef2SBarry Smith   PetscFunctionBegin;
2460089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2461089b2837SJed Brown   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
246224989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
246324989b8cSPeter Brune   ierr = DMTSGetRHSJacobian(dm,func,ctx);CHKERRQ(ierr);
2464d763cef2SBarry Smith   PetscFunctionReturn(0);
2465d763cef2SBarry Smith }
2466d763cef2SBarry Smith 
24671713a123SBarry Smith #undef __FUNCT__
24682eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian"
24692eca1d9cSJed Brown /*@C
24702eca1d9cSJed Brown    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
24712eca1d9cSJed Brown 
24722eca1d9cSJed Brown    Not Collective, but parallel objects are returned if TS is parallel
24732eca1d9cSJed Brown 
24742eca1d9cSJed Brown    Input Parameter:
24752eca1d9cSJed Brown .  ts  - The TS context obtained from TSCreate()
24762eca1d9cSJed Brown 
24772eca1d9cSJed Brown    Output Parameters:
24782eca1d9cSJed Brown +  A   - The Jacobian of F(t,U,U_t)
24792eca1d9cSJed Brown .  B   - The preconditioner matrix, often the same as A
24802eca1d9cSJed Brown .  f   - The function to compute the matrices
24812eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine
24822eca1d9cSJed Brown 
24832eca1d9cSJed Brown    Notes: You can pass in PETSC_NULL for any return argument you do not need.
24842eca1d9cSJed Brown 
24852eca1d9cSJed Brown    Level: advanced
24862eca1d9cSJed Brown 
24872eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
24882eca1d9cSJed Brown 
24892eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian
24902eca1d9cSJed Brown @*/
24917087cfbeSBarry Smith PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
24922eca1d9cSJed Brown {
2493089b2837SJed Brown   PetscErrorCode ierr;
2494089b2837SJed Brown   SNES           snes;
249524989b8cSPeter Brune   DM             dm;
24962eca1d9cSJed Brown   PetscFunctionBegin;
2497089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2498089b2837SJed Brown   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
249924989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
250024989b8cSPeter Brune   ierr = DMTSGetIJacobian(dm,f,ctx);CHKERRQ(ierr);
25012eca1d9cSJed Brown   PetscFunctionReturn(0);
25022eca1d9cSJed Brown }
25032eca1d9cSJed Brown 
25046083293cSBarry Smith typedef struct {
25056083293cSBarry Smith   PetscViewer viewer;
25066083293cSBarry Smith   Vec         initialsolution;
25076083293cSBarry Smith   PetscBool   showinitial;
25086083293cSBarry Smith } TSMonitorSolutionCtx;
25096083293cSBarry Smith 
25102eca1d9cSJed Brown #undef __FUNCT__
2511a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution"
25121713a123SBarry Smith /*@C
2513a6570f20SBarry Smith    TSMonitorSolution - Monitors progress of the TS solvers by calling
25141713a123SBarry Smith    VecView() for the solution at each timestep
25151713a123SBarry Smith 
25161713a123SBarry Smith    Collective on TS
25171713a123SBarry Smith 
25181713a123SBarry Smith    Input Parameters:
25191713a123SBarry Smith +  ts - the TS context
25201713a123SBarry Smith .  step - current time-step
2521142b95e3SSatish Balay .  ptime - current time
25221713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
25231713a123SBarry Smith 
25241713a123SBarry Smith    Level: intermediate
25251713a123SBarry Smith 
25261713a123SBarry Smith .keywords: TS,  vector, monitor, view
25271713a123SBarry Smith 
2528a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
25291713a123SBarry Smith @*/
25307087cfbeSBarry Smith PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
25311713a123SBarry Smith {
2532dfbe8321SBarry Smith   PetscErrorCode       ierr;
25336083293cSBarry Smith   TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy;
25341713a123SBarry Smith 
25351713a123SBarry Smith   PetscFunctionBegin;
25366083293cSBarry Smith   if (!step && ictx->showinitial) {
25376083293cSBarry Smith     if (!ictx->initialsolution) {
25386083293cSBarry Smith       ierr = VecDuplicate(x,&ictx->initialsolution);CHKERRQ(ierr);
25391713a123SBarry Smith     }
25406083293cSBarry Smith     ierr = VecCopy(x,ictx->initialsolution);CHKERRQ(ierr);
25416083293cSBarry Smith   }
25426083293cSBarry Smith   if (ictx->showinitial) {
25436083293cSBarry Smith     PetscReal pause;
25446083293cSBarry Smith     ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr);
25456083293cSBarry Smith     ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr);
25466083293cSBarry Smith     ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr);
25476083293cSBarry Smith     ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr);
25486083293cSBarry Smith     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr);
25496083293cSBarry Smith   }
25506083293cSBarry Smith   ierr = VecView(x,ictx->viewer);CHKERRQ(ierr);
25516083293cSBarry Smith   if (ictx->showinitial) {
25526083293cSBarry Smith     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr);
25536083293cSBarry Smith   }
25541713a123SBarry Smith   PetscFunctionReturn(0);
25551713a123SBarry Smith }
25561713a123SBarry Smith 
25571713a123SBarry Smith 
25586c699258SBarry Smith #undef __FUNCT__
25596083293cSBarry Smith #define __FUNCT__ "TSMonitorSolutionDestroy"
25606083293cSBarry Smith /*@C
25616083293cSBarry Smith    TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution
25626083293cSBarry Smith 
25636083293cSBarry Smith    Collective on TS
25646083293cSBarry Smith 
25656083293cSBarry Smith    Input Parameters:
25666083293cSBarry Smith .    ctx - the monitor context
25676083293cSBarry Smith 
25686083293cSBarry Smith    Level: intermediate
25696083293cSBarry Smith 
25706083293cSBarry Smith .keywords: TS,  vector, monitor, view
25716083293cSBarry Smith 
25726083293cSBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
25736083293cSBarry Smith @*/
25746083293cSBarry Smith PetscErrorCode  TSMonitorSolutionDestroy(void **ctx)
25756083293cSBarry Smith {
25766083293cSBarry Smith   PetscErrorCode       ierr;
25776083293cSBarry Smith   TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx;
25786083293cSBarry Smith 
25796083293cSBarry Smith   PetscFunctionBegin;
25806083293cSBarry Smith   ierr = PetscViewerDestroy(&ictx->viewer);CHKERRQ(ierr);
25816083293cSBarry Smith   ierr = VecDestroy(&ictx->initialsolution);CHKERRQ(ierr);
25826083293cSBarry Smith   ierr = PetscFree(ictx);CHKERRQ(ierr);
25836083293cSBarry Smith   PetscFunctionReturn(0);
25846083293cSBarry Smith }
25856083293cSBarry Smith 
25866083293cSBarry Smith #undef __FUNCT__
25876083293cSBarry Smith #define __FUNCT__ "TSMonitorSolutionCreate"
25886083293cSBarry Smith /*@C
25896083293cSBarry Smith    TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution
25906083293cSBarry Smith 
25916083293cSBarry Smith    Collective on TS
25926083293cSBarry Smith 
25936083293cSBarry Smith    Input Parameter:
25946083293cSBarry Smith .    ts - time-step context
25956083293cSBarry Smith 
25966083293cSBarry Smith    Output Patameter:
25976083293cSBarry Smith .    ctx - the monitor context
25986083293cSBarry Smith 
25996083293cSBarry Smith    Level: intermediate
26006083293cSBarry Smith 
26016083293cSBarry Smith .keywords: TS,  vector, monitor, view
26026083293cSBarry Smith 
26036083293cSBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
26046083293cSBarry Smith @*/
26056083293cSBarry Smith PetscErrorCode  TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx)
26066083293cSBarry Smith {
26076083293cSBarry Smith   PetscErrorCode       ierr;
26086083293cSBarry Smith   TSMonitorSolutionCtx *ictx;
26096083293cSBarry Smith 
26106083293cSBarry Smith   PetscFunctionBegin;
26116083293cSBarry Smith   ierr = PetscNew(TSMonitorSolutionCtx,&ictx);CHKERRQ(ierr);
26126083293cSBarry Smith   *ctx = (void*)ictx;
26136083293cSBarry Smith   if (!viewer) {
26146083293cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
26156083293cSBarry Smith   }
26166083293cSBarry Smith   ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
26176083293cSBarry Smith   ictx->viewer      = viewer;
26186083293cSBarry Smith   ictx->showinitial = PETSC_FALSE;
26196083293cSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);CHKERRQ(ierr);
26206083293cSBarry Smith   PetscFunctionReturn(0);
26216083293cSBarry Smith }
26226083293cSBarry Smith 
26236083293cSBarry Smith #undef __FUNCT__
26246c699258SBarry Smith #define __FUNCT__ "TSSetDM"
26256c699258SBarry Smith /*@
26266c699258SBarry Smith    TSSetDM - Sets the DM that may be used by some preconditioners
26276c699258SBarry Smith 
26283f9fe445SBarry Smith    Logically Collective on TS and DM
26296c699258SBarry Smith 
26306c699258SBarry Smith    Input Parameters:
26316c699258SBarry Smith +  ts - the preconditioner context
26326c699258SBarry Smith -  dm - the dm
26336c699258SBarry Smith 
26346c699258SBarry Smith    Level: intermediate
26356c699258SBarry Smith 
26366c699258SBarry Smith 
26376c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
26386c699258SBarry Smith @*/
26397087cfbeSBarry Smith PetscErrorCode  TSSetDM(TS ts,DM dm)
26406c699258SBarry Smith {
26416c699258SBarry Smith   PetscErrorCode ierr;
2642089b2837SJed Brown   SNES           snes;
264324989b8cSPeter Brune   TSDM           tsdm;
26446c699258SBarry Smith 
26456c699258SBarry Smith   PetscFunctionBegin;
26460700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
264770663e4aSLisandro Dalcin   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
264824989b8cSPeter Brune   if (ts->dm) {               /* Move the TSDM context over to the new DM unless the new DM already has one */
264924989b8cSPeter Brune     PetscContainer oldcontainer,container;
265024989b8cSPeter Brune     ierr = PetscObjectQuery((PetscObject)ts->dm,"TSDM",(PetscObject*)&oldcontainer);CHKERRQ(ierr);
265124989b8cSPeter Brune     ierr = PetscObjectQuery((PetscObject)dm,"TSDM",(PetscObject*)&container);CHKERRQ(ierr);
265224989b8cSPeter Brune     if (oldcontainer && !container) {
265324989b8cSPeter Brune       ierr = DMTSCopyContext(ts->dm,dm);CHKERRQ(ierr);
265424989b8cSPeter Brune       ierr = DMTSGetContext(ts->dm,&tsdm);CHKERRQ(ierr);
265524989b8cSPeter Brune       if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
265624989b8cSPeter Brune         tsdm->originaldm = dm;
265724989b8cSPeter Brune       }
265824989b8cSPeter Brune     }
26596bf464f9SBarry Smith     ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
266024989b8cSPeter Brune   }
26616c699258SBarry Smith   ts->dm = dm;
2662089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2663089b2837SJed Brown   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
26646c699258SBarry Smith   PetscFunctionReturn(0);
26656c699258SBarry Smith }
26666c699258SBarry Smith 
26676c699258SBarry Smith #undef __FUNCT__
26686c699258SBarry Smith #define __FUNCT__ "TSGetDM"
26696c699258SBarry Smith /*@
26706c699258SBarry Smith    TSGetDM - Gets the DM that may be used by some preconditioners
26716c699258SBarry Smith 
26723f9fe445SBarry Smith    Not Collective
26736c699258SBarry Smith 
26746c699258SBarry Smith    Input Parameter:
26756c699258SBarry Smith . ts - the preconditioner context
26766c699258SBarry Smith 
26776c699258SBarry Smith    Output Parameter:
26786c699258SBarry Smith .  dm - the dm
26796c699258SBarry Smith 
26806c699258SBarry Smith    Level: intermediate
26816c699258SBarry Smith 
26826c699258SBarry Smith 
26836c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
26846c699258SBarry Smith @*/
26857087cfbeSBarry Smith PetscErrorCode  TSGetDM(TS ts,DM *dm)
26866c699258SBarry Smith {
2687496e6a7aSJed Brown   PetscErrorCode ierr;
2688496e6a7aSJed Brown 
26896c699258SBarry Smith   PetscFunctionBegin;
26900700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2691496e6a7aSJed Brown   if (!ts->dm) {
2692496e6a7aSJed Brown     ierr = DMShellCreate(((PetscObject)ts)->comm,&ts->dm);CHKERRQ(ierr);
2693496e6a7aSJed Brown     if (ts->snes) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
2694496e6a7aSJed Brown   }
26956c699258SBarry Smith   *dm = ts->dm;
26966c699258SBarry Smith   PetscFunctionReturn(0);
26976c699258SBarry Smith }
26981713a123SBarry Smith 
26990f5c6efeSJed Brown #undef __FUNCT__
27000f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction"
27010f5c6efeSJed Brown /*@
27020f5c6efeSJed Brown    SNESTSFormFunction - Function to evaluate nonlinear residual
27030f5c6efeSJed Brown 
27043f9fe445SBarry Smith    Logically Collective on SNES
27050f5c6efeSJed Brown 
27060f5c6efeSJed Brown    Input Parameter:
2707d42a1c89SJed Brown + snes - nonlinear solver
27080f5c6efeSJed Brown . X - the current state at which to evaluate the residual
2709d42a1c89SJed Brown - ctx - user context, must be a TS
27100f5c6efeSJed Brown 
27110f5c6efeSJed Brown    Output Parameter:
27120f5c6efeSJed Brown . F - the nonlinear residual
27130f5c6efeSJed Brown 
27140f5c6efeSJed Brown    Notes:
27150f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
27160f5c6efeSJed Brown    It is most frequently passed to MatFDColoringSetFunction().
27170f5c6efeSJed Brown 
27180f5c6efeSJed Brown    Level: advanced
27190f5c6efeSJed Brown 
27200f5c6efeSJed Brown .seealso: SNESSetFunction(), MatFDColoringSetFunction()
27210f5c6efeSJed Brown @*/
27227087cfbeSBarry Smith PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
27230f5c6efeSJed Brown {
27240f5c6efeSJed Brown   TS ts = (TS)ctx;
27250f5c6efeSJed Brown   PetscErrorCode ierr;
27260f5c6efeSJed Brown 
27270f5c6efeSJed Brown   PetscFunctionBegin;
27280f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
27290f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
27300f5c6efeSJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
27310f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
27320f5c6efeSJed Brown   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
27330f5c6efeSJed Brown   PetscFunctionReturn(0);
27340f5c6efeSJed Brown }
27350f5c6efeSJed Brown 
27360f5c6efeSJed Brown #undef __FUNCT__
27370f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian"
27380f5c6efeSJed Brown /*@
27390f5c6efeSJed Brown    SNESTSFormJacobian - Function to evaluate the Jacobian
27400f5c6efeSJed Brown 
27410f5c6efeSJed Brown    Collective on SNES
27420f5c6efeSJed Brown 
27430f5c6efeSJed Brown    Input Parameter:
27440f5c6efeSJed Brown + snes - nonlinear solver
27450f5c6efeSJed Brown . X - the current state at which to evaluate the residual
27460f5c6efeSJed Brown - ctx - user context, must be a TS
27470f5c6efeSJed Brown 
27480f5c6efeSJed Brown    Output Parameter:
27490f5c6efeSJed Brown + A - the Jacobian
27500f5c6efeSJed Brown . B - the preconditioning matrix (may be the same as A)
27510f5c6efeSJed Brown - flag - indicates any structure change in the matrix
27520f5c6efeSJed Brown 
27530f5c6efeSJed Brown    Notes:
27540f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
27550f5c6efeSJed Brown 
27560f5c6efeSJed Brown    Level: developer
27570f5c6efeSJed Brown 
27580f5c6efeSJed Brown .seealso: SNESSetJacobian()
27590f5c6efeSJed Brown @*/
27607087cfbeSBarry Smith PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
27610f5c6efeSJed Brown {
27620f5c6efeSJed Brown   TS ts = (TS)ctx;
27630f5c6efeSJed Brown   PetscErrorCode ierr;
27640f5c6efeSJed Brown 
27650f5c6efeSJed Brown   PetscFunctionBegin;
27660f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
27670f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
27680f5c6efeSJed Brown   PetscValidPointer(A,3);
27690f5c6efeSJed Brown   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
27700f5c6efeSJed Brown   PetscValidPointer(B,4);
27710f5c6efeSJed Brown   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
27720f5c6efeSJed Brown   PetscValidPointer(flag,5);
27730f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
27740f5c6efeSJed Brown   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
27750f5c6efeSJed Brown   PetscFunctionReturn(0);
27760f5c6efeSJed Brown }
2777325fc9f4SBarry Smith 
27780e4ef248SJed Brown #undef __FUNCT__
27790e4ef248SJed Brown #define __FUNCT__ "TSComputeRHSFunctionLinear"
27800e4ef248SJed Brown /*@C
27810e4ef248SJed Brown    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only
27820e4ef248SJed Brown 
27830e4ef248SJed Brown    Collective on TS
27840e4ef248SJed Brown 
27850e4ef248SJed Brown    Input Arguments:
27860e4ef248SJed Brown +  ts - time stepping context
27870e4ef248SJed Brown .  t - time at which to evaluate
27880e4ef248SJed Brown .  X - state at which to evaluate
27890e4ef248SJed Brown -  ctx - context
27900e4ef248SJed Brown 
27910e4ef248SJed Brown    Output Arguments:
27920e4ef248SJed Brown .  F - right hand side
27930e4ef248SJed Brown 
27940e4ef248SJed Brown    Level: intermediate
27950e4ef248SJed Brown 
27960e4ef248SJed Brown    Notes:
27970e4ef248SJed Brown    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
27980e4ef248SJed Brown    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
27990e4ef248SJed Brown 
28000e4ef248SJed Brown .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
28010e4ef248SJed Brown @*/
28020e4ef248SJed Brown PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
28030e4ef248SJed Brown {
28040e4ef248SJed Brown   PetscErrorCode ierr;
28050e4ef248SJed Brown   Mat Arhs,Brhs;
28060e4ef248SJed Brown   MatStructure flg2;
28070e4ef248SJed Brown 
28080e4ef248SJed Brown   PetscFunctionBegin;
28090e4ef248SJed Brown   ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
28100e4ef248SJed Brown   ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
28110e4ef248SJed Brown   ierr = MatMult(Arhs,X,F);CHKERRQ(ierr);
28120e4ef248SJed Brown   PetscFunctionReturn(0);
28130e4ef248SJed Brown }
28140e4ef248SJed Brown 
28150e4ef248SJed Brown #undef __FUNCT__
28160e4ef248SJed Brown #define __FUNCT__ "TSComputeRHSJacobianConstant"
28170e4ef248SJed Brown /*@C
28180e4ef248SJed Brown    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
28190e4ef248SJed Brown 
28200e4ef248SJed Brown    Collective on TS
28210e4ef248SJed Brown 
28220e4ef248SJed Brown    Input Arguments:
28230e4ef248SJed Brown +  ts - time stepping context
28240e4ef248SJed Brown .  t - time at which to evaluate
28250e4ef248SJed Brown .  X - state at which to evaluate
28260e4ef248SJed Brown -  ctx - context
28270e4ef248SJed Brown 
28280e4ef248SJed Brown    Output Arguments:
28290e4ef248SJed Brown +  A - pointer to operator
28300e4ef248SJed Brown .  B - pointer to preconditioning matrix
28310e4ef248SJed Brown -  flg - matrix structure flag
28320e4ef248SJed Brown 
28330e4ef248SJed Brown    Level: intermediate
28340e4ef248SJed Brown 
28350e4ef248SJed Brown    Notes:
28360e4ef248SJed Brown    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
28370e4ef248SJed Brown 
28380e4ef248SJed Brown .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
28390e4ef248SJed Brown @*/
28400e4ef248SJed Brown PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
28410e4ef248SJed Brown {
28420e4ef248SJed Brown 
28430e4ef248SJed Brown   PetscFunctionBegin;
28440e4ef248SJed Brown   *flg = SAME_PRECONDITIONER;
28450e4ef248SJed Brown   PetscFunctionReturn(0);
28460e4ef248SJed Brown }
28470e4ef248SJed Brown 
28480026cea9SSean Farley #undef __FUNCT__
28490026cea9SSean Farley #define __FUNCT__ "TSComputeIFunctionLinear"
28500026cea9SSean Farley /*@C
28510026cea9SSean Farley    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
28520026cea9SSean Farley 
28530026cea9SSean Farley    Collective on TS
28540026cea9SSean Farley 
28550026cea9SSean Farley    Input Arguments:
28560026cea9SSean Farley +  ts - time stepping context
28570026cea9SSean Farley .  t - time at which to evaluate
28580026cea9SSean Farley .  X - state at which to evaluate
28590026cea9SSean Farley .  Xdot - time derivative of state vector
28600026cea9SSean Farley -  ctx - context
28610026cea9SSean Farley 
28620026cea9SSean Farley    Output Arguments:
28630026cea9SSean Farley .  F - left hand side
28640026cea9SSean Farley 
28650026cea9SSean Farley    Level: intermediate
28660026cea9SSean Farley 
28670026cea9SSean Farley    Notes:
28680026cea9SSean 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
28690026cea9SSean Farley    user is required to write their own TSComputeIFunction.
28700026cea9SSean Farley    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
28710026cea9SSean Farley    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
28720026cea9SSean Farley 
28730026cea9SSean Farley .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
28740026cea9SSean Farley @*/
28750026cea9SSean Farley PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
28760026cea9SSean Farley {
28770026cea9SSean Farley   PetscErrorCode ierr;
28780026cea9SSean Farley   Mat A,B;
28790026cea9SSean Farley   MatStructure flg2;
28800026cea9SSean Farley 
28810026cea9SSean Farley   PetscFunctionBegin;
28820026cea9SSean Farley   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
28830026cea9SSean Farley   ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr);
28840026cea9SSean Farley   ierr = MatMult(A,Xdot,F);CHKERRQ(ierr);
28850026cea9SSean Farley   PetscFunctionReturn(0);
28860026cea9SSean Farley }
28870026cea9SSean Farley 
28880026cea9SSean Farley #undef __FUNCT__
28890026cea9SSean Farley #define __FUNCT__ "TSComputeIJacobianConstant"
28900026cea9SSean Farley /*@C
289124e94211SJed Brown    TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent.
28920026cea9SSean Farley 
28930026cea9SSean Farley    Collective on TS
28940026cea9SSean Farley 
28950026cea9SSean Farley    Input Arguments:
28960026cea9SSean Farley +  ts - time stepping context
28970026cea9SSean Farley .  t - time at which to evaluate
28980026cea9SSean Farley .  X - state at which to evaluate
28990026cea9SSean Farley .  Xdot - time derivative of state vector
29000026cea9SSean Farley .  shift - shift to apply
29010026cea9SSean Farley -  ctx - context
29020026cea9SSean Farley 
29030026cea9SSean Farley    Output Arguments:
29040026cea9SSean Farley +  A - pointer to operator
29050026cea9SSean Farley .  B - pointer to preconditioning matrix
29060026cea9SSean Farley -  flg - matrix structure flag
29070026cea9SSean Farley 
29080026cea9SSean Farley    Level: intermediate
29090026cea9SSean Farley 
29100026cea9SSean Farley    Notes:
29110026cea9SSean Farley    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
29120026cea9SSean Farley 
29130026cea9SSean Farley .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
29140026cea9SSean Farley @*/
29150026cea9SSean Farley PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
29160026cea9SSean Farley {
29170026cea9SSean Farley 
29180026cea9SSean Farley   PetscFunctionBegin;
29190026cea9SSean Farley   *flg = SAME_PRECONDITIONER;
29200026cea9SSean Farley   PetscFunctionReturn(0);
29210026cea9SSean Farley }
29220026cea9SSean Farley 
29234af1b03aSJed Brown 
29244af1b03aSJed Brown #undef __FUNCT__
29254af1b03aSJed Brown #define __FUNCT__ "TSGetConvergedReason"
29264af1b03aSJed Brown /*@
29274af1b03aSJed Brown    TSGetConvergedReason - Gets the reason the TS iteration was stopped.
29284af1b03aSJed Brown 
29294af1b03aSJed Brown    Not Collective
29304af1b03aSJed Brown 
29314af1b03aSJed Brown    Input Parameter:
29324af1b03aSJed Brown .  ts - the TS context
29334af1b03aSJed Brown 
29344af1b03aSJed Brown    Output Parameter:
29354af1b03aSJed Brown .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
29364af1b03aSJed Brown             manual pages for the individual convergence tests for complete lists
29374af1b03aSJed Brown 
29384af1b03aSJed Brown    Level: intermediate
29394af1b03aSJed Brown 
2940cd652676SJed Brown    Notes:
2941cd652676SJed Brown    Can only be called after the call to TSSolve() is complete.
29424af1b03aSJed Brown 
29434af1b03aSJed Brown .keywords: TS, nonlinear, set, convergence, test
29444af1b03aSJed Brown 
29454af1b03aSJed Brown .seealso: TSSetConvergenceTest(), TSConvergedReason
29464af1b03aSJed Brown @*/
29474af1b03aSJed Brown PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
29484af1b03aSJed Brown {
29494af1b03aSJed Brown   PetscFunctionBegin;
29504af1b03aSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
29514af1b03aSJed Brown   PetscValidPointer(reason,2);
29524af1b03aSJed Brown   *reason = ts->reason;
29534af1b03aSJed Brown   PetscFunctionReturn(0);
29544af1b03aSJed Brown }
29554af1b03aSJed Brown 
2956fb1732b5SBarry Smith #undef __FUNCT__
29575ef26d82SJed Brown #define __FUNCT__ "TSGetSNESIterations"
29589f67acb7SJed Brown /*@
29595ef26d82SJed Brown    TSGetSNESIterations - Gets the total number of nonlinear iterations
29609f67acb7SJed Brown    used by the time integrator.
29619f67acb7SJed Brown 
29629f67acb7SJed Brown    Not Collective
29639f67acb7SJed Brown 
29649f67acb7SJed Brown    Input Parameter:
29659f67acb7SJed Brown .  ts - TS context
29669f67acb7SJed Brown 
29679f67acb7SJed Brown    Output Parameter:
29689f67acb7SJed Brown .  nits - number of nonlinear iterations
29699f67acb7SJed Brown 
29709f67acb7SJed Brown    Notes:
29719f67acb7SJed Brown    This counter is reset to zero for each successive call to TSSolve().
29729f67acb7SJed Brown 
29739f67acb7SJed Brown    Level: intermediate
29749f67acb7SJed Brown 
29759f67acb7SJed Brown .keywords: TS, get, number, nonlinear, iterations
29769f67acb7SJed Brown 
29775ef26d82SJed Brown .seealso:  TSGetKSPIterations()
29789f67acb7SJed Brown @*/
29795ef26d82SJed Brown PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
29809f67acb7SJed Brown {
29819f67acb7SJed Brown   PetscFunctionBegin;
29829f67acb7SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
29839f67acb7SJed Brown   PetscValidIntPointer(nits,2);
29845ef26d82SJed Brown   *nits = ts->snes_its;
29859f67acb7SJed Brown   PetscFunctionReturn(0);
29869f67acb7SJed Brown }
29879f67acb7SJed Brown 
29889f67acb7SJed Brown #undef __FUNCT__
29895ef26d82SJed Brown #define __FUNCT__ "TSGetKSPIterations"
29909f67acb7SJed Brown /*@
29915ef26d82SJed Brown    TSGetKSPIterations - Gets the total number of linear iterations
29929f67acb7SJed Brown    used by the time integrator.
29939f67acb7SJed Brown 
29949f67acb7SJed Brown    Not Collective
29959f67acb7SJed Brown 
29969f67acb7SJed Brown    Input Parameter:
29979f67acb7SJed Brown .  ts - TS context
29989f67acb7SJed Brown 
29999f67acb7SJed Brown    Output Parameter:
30009f67acb7SJed Brown .  lits - number of linear iterations
30019f67acb7SJed Brown 
30029f67acb7SJed Brown    Notes:
30039f67acb7SJed Brown    This counter is reset to zero for each successive call to TSSolve().
30049f67acb7SJed Brown 
30059f67acb7SJed Brown    Level: intermediate
30069f67acb7SJed Brown 
30079f67acb7SJed Brown .keywords: TS, get, number, linear, iterations
30089f67acb7SJed Brown 
30095ef26d82SJed Brown .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
30109f67acb7SJed Brown @*/
30115ef26d82SJed Brown PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
30129f67acb7SJed Brown {
30139f67acb7SJed Brown   PetscFunctionBegin;
30149f67acb7SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
30159f67acb7SJed Brown   PetscValidIntPointer(lits,2);
30165ef26d82SJed Brown   *lits = ts->ksp_its;
30179f67acb7SJed Brown   PetscFunctionReturn(0);
30189f67acb7SJed Brown }
30199f67acb7SJed Brown 
30209f67acb7SJed Brown #undef __FUNCT__
3021cef5090cSJed Brown #define __FUNCT__ "TSGetStepRejections"
3022cef5090cSJed Brown /*@
3023cef5090cSJed Brown    TSGetStepRejections - Gets the total number of rejected steps.
3024cef5090cSJed Brown 
3025cef5090cSJed Brown    Not Collective
3026cef5090cSJed Brown 
3027cef5090cSJed Brown    Input Parameter:
3028cef5090cSJed Brown .  ts - TS context
3029cef5090cSJed Brown 
3030cef5090cSJed Brown    Output Parameter:
3031cef5090cSJed Brown .  rejects - number of steps rejected
3032cef5090cSJed Brown 
3033cef5090cSJed Brown    Notes:
3034cef5090cSJed Brown    This counter is reset to zero for each successive call to TSSolve().
3035cef5090cSJed Brown 
3036cef5090cSJed Brown    Level: intermediate
3037cef5090cSJed Brown 
3038cef5090cSJed Brown .keywords: TS, get, number
3039cef5090cSJed Brown 
30405ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
3041cef5090cSJed Brown @*/
3042cef5090cSJed Brown PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
3043cef5090cSJed Brown {
3044cef5090cSJed Brown   PetscFunctionBegin;
3045cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3046cef5090cSJed Brown   PetscValidIntPointer(rejects,2);
3047cef5090cSJed Brown   *rejects = ts->reject;
3048cef5090cSJed Brown   PetscFunctionReturn(0);
3049cef5090cSJed Brown }
3050cef5090cSJed Brown 
3051cef5090cSJed Brown #undef __FUNCT__
3052cef5090cSJed Brown #define __FUNCT__ "TSGetSNESFailures"
3053cef5090cSJed Brown /*@
3054cef5090cSJed Brown    TSGetSNESFailures - Gets the total number of failed SNES solves
3055cef5090cSJed Brown 
3056cef5090cSJed Brown    Not Collective
3057cef5090cSJed Brown 
3058cef5090cSJed Brown    Input Parameter:
3059cef5090cSJed Brown .  ts - TS context
3060cef5090cSJed Brown 
3061cef5090cSJed Brown    Output Parameter:
3062cef5090cSJed Brown .  fails - number of failed nonlinear solves
3063cef5090cSJed Brown 
3064cef5090cSJed Brown    Notes:
3065cef5090cSJed Brown    This counter is reset to zero for each successive call to TSSolve().
3066cef5090cSJed Brown 
3067cef5090cSJed Brown    Level: intermediate
3068cef5090cSJed Brown 
3069cef5090cSJed Brown .keywords: TS, get, number
3070cef5090cSJed Brown 
30715ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
3072cef5090cSJed Brown @*/
3073cef5090cSJed Brown PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
3074cef5090cSJed Brown {
3075cef5090cSJed Brown   PetscFunctionBegin;
3076cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3077cef5090cSJed Brown   PetscValidIntPointer(fails,2);
3078cef5090cSJed Brown   *fails = ts->num_snes_failures;
3079cef5090cSJed Brown   PetscFunctionReturn(0);
3080cef5090cSJed Brown }
3081cef5090cSJed Brown 
3082cef5090cSJed Brown #undef __FUNCT__
3083cef5090cSJed Brown #define __FUNCT__ "TSSetMaxStepRejections"
3084cef5090cSJed Brown /*@
3085cef5090cSJed Brown    TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails
3086cef5090cSJed Brown 
3087cef5090cSJed Brown    Not Collective
3088cef5090cSJed Brown 
3089cef5090cSJed Brown    Input Parameter:
3090cef5090cSJed Brown +  ts - TS context
3091cef5090cSJed Brown -  rejects - maximum number of rejected steps, pass -1 for unlimited
3092cef5090cSJed Brown 
3093cef5090cSJed Brown    Notes:
3094cef5090cSJed Brown    The counter is reset to zero for each step
3095cef5090cSJed Brown 
3096cef5090cSJed Brown    Options Database Key:
3097cef5090cSJed Brown  .  -ts_max_reject - Maximum number of step rejections before a step fails
3098cef5090cSJed Brown 
3099cef5090cSJed Brown    Level: intermediate
3100cef5090cSJed Brown 
3101cef5090cSJed Brown .keywords: TS, set, maximum, number
3102cef5090cSJed Brown 
31035ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3104cef5090cSJed Brown @*/
3105cef5090cSJed Brown PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
3106cef5090cSJed Brown {
3107cef5090cSJed Brown   PetscFunctionBegin;
3108cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3109cef5090cSJed Brown   ts->max_reject = rejects;
3110cef5090cSJed Brown   PetscFunctionReturn(0);
3111cef5090cSJed Brown }
3112cef5090cSJed Brown 
3113cef5090cSJed Brown #undef __FUNCT__
3114cef5090cSJed Brown #define __FUNCT__ "TSSetMaxSNESFailures"
3115cef5090cSJed Brown /*@
3116cef5090cSJed Brown    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves
3117cef5090cSJed Brown 
3118cef5090cSJed Brown    Not Collective
3119cef5090cSJed Brown 
3120cef5090cSJed Brown    Input Parameter:
3121cef5090cSJed Brown +  ts - TS context
3122cef5090cSJed Brown -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited
3123cef5090cSJed Brown 
3124cef5090cSJed Brown    Notes:
3125cef5090cSJed Brown    The counter is reset to zero for each successive call to TSSolve().
3126cef5090cSJed Brown 
3127cef5090cSJed Brown    Options Database Key:
3128cef5090cSJed Brown  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures
3129cef5090cSJed Brown 
3130cef5090cSJed Brown    Level: intermediate
3131cef5090cSJed Brown 
3132cef5090cSJed Brown .keywords: TS, set, maximum, number
3133cef5090cSJed Brown 
31345ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
3135cef5090cSJed Brown @*/
3136cef5090cSJed Brown PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
3137cef5090cSJed Brown {
3138cef5090cSJed Brown   PetscFunctionBegin;
3139cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3140cef5090cSJed Brown   ts->max_snes_failures = fails;
3141cef5090cSJed Brown   PetscFunctionReturn(0);
3142cef5090cSJed Brown }
3143cef5090cSJed Brown 
3144cef5090cSJed Brown #undef __FUNCT__
3145cef5090cSJed Brown #define __FUNCT__ "TSSetErrorIfStepFails()"
3146cef5090cSJed Brown /*@
3147cef5090cSJed Brown    TSSetErrorIfStepFails - Error if no step succeeds
3148cef5090cSJed Brown 
3149cef5090cSJed Brown    Not Collective
3150cef5090cSJed Brown 
3151cef5090cSJed Brown    Input Parameter:
3152cef5090cSJed Brown +  ts - TS context
3153cef5090cSJed Brown -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure
3154cef5090cSJed Brown 
3155cef5090cSJed Brown    Options Database Key:
3156cef5090cSJed Brown  .  -ts_error_if_step_fails - Error if no step succeeds
3157cef5090cSJed Brown 
3158cef5090cSJed Brown    Level: intermediate
3159cef5090cSJed Brown 
3160cef5090cSJed Brown .keywords: TS, set, error
3161cef5090cSJed Brown 
31625ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3163cef5090cSJed Brown @*/
3164cef5090cSJed Brown PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
3165cef5090cSJed Brown {
3166cef5090cSJed Brown   PetscFunctionBegin;
3167cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3168cef5090cSJed Brown   ts->errorifstepfailed = err;
3169cef5090cSJed Brown   PetscFunctionReturn(0);
3170cef5090cSJed Brown }
3171cef5090cSJed Brown 
3172cef5090cSJed Brown #undef __FUNCT__
3173fb1732b5SBarry Smith #define __FUNCT__ "TSMonitorSolutionBinary"
3174fb1732b5SBarry Smith /*@C
3175fb1732b5SBarry Smith    TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file
3176fb1732b5SBarry Smith 
3177fb1732b5SBarry Smith    Collective on TS
3178fb1732b5SBarry Smith 
3179fb1732b5SBarry Smith    Input Parameters:
3180fb1732b5SBarry Smith +  ts - the TS context
3181fb1732b5SBarry Smith .  step - current time-step
3182fb1732b5SBarry Smith .  ptime - current time
3183ed81e22dSJed Brown .  x - current state
3184fb1732b5SBarry Smith -  viewer - binary viewer
3185fb1732b5SBarry Smith 
3186fb1732b5SBarry Smith    Level: intermediate
3187fb1732b5SBarry Smith 
3188fb1732b5SBarry Smith .keywords: TS,  vector, monitor, view
3189fb1732b5SBarry Smith 
3190fb1732b5SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3191fb1732b5SBarry Smith @*/
3192ed81e22dSJed Brown PetscErrorCode  TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer)
3193fb1732b5SBarry Smith {
3194fb1732b5SBarry Smith   PetscErrorCode       ierr;
3195ed81e22dSJed Brown   PetscViewer          v = (PetscViewer)viewer;
3196fb1732b5SBarry Smith 
3197fb1732b5SBarry Smith   PetscFunctionBegin;
3198ed81e22dSJed Brown   ierr = VecView(x,v);CHKERRQ(ierr);
3199ed81e22dSJed Brown   PetscFunctionReturn(0);
3200ed81e22dSJed Brown }
3201ed81e22dSJed Brown 
3202ed81e22dSJed Brown #undef __FUNCT__
3203ed81e22dSJed Brown #define __FUNCT__ "TSMonitorSolutionVTK"
3204ed81e22dSJed Brown /*@C
3205ed81e22dSJed Brown    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.
3206ed81e22dSJed Brown 
3207ed81e22dSJed Brown    Collective on TS
3208ed81e22dSJed Brown 
3209ed81e22dSJed Brown    Input Parameters:
3210ed81e22dSJed Brown +  ts - the TS context
3211ed81e22dSJed Brown .  step - current time-step
3212ed81e22dSJed Brown .  ptime - current time
3213ed81e22dSJed Brown .  x - current state
3214ed81e22dSJed Brown -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3215ed81e22dSJed Brown 
3216ed81e22dSJed Brown    Level: intermediate
3217ed81e22dSJed Brown 
3218ed81e22dSJed Brown    Notes:
3219ed81e22dSJed Brown    The VTK format does not allow writing multiple time steps in the same file, therefore a different file will be written for each time step.
3220ed81e22dSJed Brown    These are named according to the file name template.
3221ed81e22dSJed Brown 
3222ed81e22dSJed Brown    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().
3223ed81e22dSJed Brown 
3224ed81e22dSJed Brown .keywords: TS,  vector, monitor, view
3225ed81e22dSJed Brown 
3226ed81e22dSJed Brown .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3227ed81e22dSJed Brown @*/
3228ed81e22dSJed Brown PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate)
3229ed81e22dSJed Brown {
3230ed81e22dSJed Brown   PetscErrorCode ierr;
3231ed81e22dSJed Brown   char           filename[PETSC_MAX_PATH_LEN];
3232ed81e22dSJed Brown   PetscViewer    viewer;
3233ed81e22dSJed Brown 
3234ed81e22dSJed Brown   PetscFunctionBegin;
3235*8caf3d72SBarry Smith   ierr = PetscSNPrintf(filename,sizeof(filename),(const char*)filenametemplate,step);CHKERRQ(ierr);
3236ed81e22dSJed Brown   ierr = PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
3237fb1732b5SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
3238ed81e22dSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3239ed81e22dSJed Brown   PetscFunctionReturn(0);
3240ed81e22dSJed Brown }
3241ed81e22dSJed Brown 
3242ed81e22dSJed Brown #undef __FUNCT__
3243ed81e22dSJed Brown #define __FUNCT__ "TSMonitorSolutionVTKDestroy"
3244ed81e22dSJed Brown /*@C
3245ed81e22dSJed Brown    TSMonitorSolutionVTKDestroy - Destroy context for monitoring
3246ed81e22dSJed Brown 
3247ed81e22dSJed Brown    Collective on TS
3248ed81e22dSJed Brown 
3249ed81e22dSJed Brown    Input Parameters:
3250ed81e22dSJed Brown .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3251ed81e22dSJed Brown 
3252ed81e22dSJed Brown    Level: intermediate
3253ed81e22dSJed Brown 
3254ed81e22dSJed Brown    Note:
3255ed81e22dSJed Brown    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().
3256ed81e22dSJed Brown 
3257ed81e22dSJed Brown .keywords: TS,  vector, monitor, view
3258ed81e22dSJed Brown 
3259ed81e22dSJed Brown .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
3260ed81e22dSJed Brown @*/
3261ed81e22dSJed Brown PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
3262ed81e22dSJed Brown {
3263ed81e22dSJed Brown   PetscErrorCode ierr;
3264ed81e22dSJed Brown 
3265ed81e22dSJed Brown   PetscFunctionBegin;
3266ed81e22dSJed Brown   ierr = PetscFree(*(char**)filenametemplate);CHKERRQ(ierr);
3267fb1732b5SBarry Smith   PetscFunctionReturn(0);
3268fb1732b5SBarry Smith }
3269fb1732b5SBarry Smith 
327084df9cb4SJed Brown #undef __FUNCT__
327184df9cb4SJed Brown #define __FUNCT__ "TSGetAdapt"
327284df9cb4SJed Brown /*@
327384df9cb4SJed Brown    TSGetAdapt - Get the adaptive controller context for the current method
327484df9cb4SJed Brown 
3275ed81e22dSJed Brown    Collective on TS if controller has not been created yet
327684df9cb4SJed Brown 
327784df9cb4SJed Brown    Input Arguments:
3278ed81e22dSJed Brown .  ts - time stepping context
327984df9cb4SJed Brown 
328084df9cb4SJed Brown    Output Arguments:
3281ed81e22dSJed Brown .  adapt - adaptive controller
328284df9cb4SJed Brown 
328384df9cb4SJed Brown    Level: intermediate
328484df9cb4SJed Brown 
3285ed81e22dSJed Brown .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
328684df9cb4SJed Brown @*/
328784df9cb4SJed Brown PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
328884df9cb4SJed Brown {
328984df9cb4SJed Brown   PetscErrorCode ierr;
329084df9cb4SJed Brown 
329184df9cb4SJed Brown   PetscFunctionBegin;
329284df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
329384df9cb4SJed Brown   PetscValidPointer(adapt,2);
329484df9cb4SJed Brown   if (!ts->adapt) {
329584df9cb4SJed Brown     ierr = TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);CHKERRQ(ierr);
32961c3436cfSJed Brown     ierr = PetscLogObjectParent(ts,ts->adapt);CHKERRQ(ierr);
32971c3436cfSJed Brown     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);CHKERRQ(ierr);
329884df9cb4SJed Brown   }
329984df9cb4SJed Brown   *adapt = ts->adapt;
330084df9cb4SJed Brown   PetscFunctionReturn(0);
330184df9cb4SJed Brown }
3302d6ebe24aSShri Abhyankar 
3303d6ebe24aSShri Abhyankar #undef __FUNCT__
33041c3436cfSJed Brown #define __FUNCT__ "TSSetTolerances"
33051c3436cfSJed Brown /*@
33061c3436cfSJed Brown    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller
33071c3436cfSJed Brown 
33081c3436cfSJed Brown    Logically Collective
33091c3436cfSJed Brown 
33101c3436cfSJed Brown    Input Arguments:
33111c3436cfSJed Brown +  ts - time integration context
33121c3436cfSJed Brown .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
33131c3436cfSJed Brown .  vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present
33141c3436cfSJed Brown .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
33151c3436cfSJed Brown -  vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present
33161c3436cfSJed Brown 
33171c3436cfSJed Brown    Level: beginner
33181c3436cfSJed Brown 
3319c5033834SJed Brown .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
33201c3436cfSJed Brown @*/
33211c3436cfSJed Brown PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
33221c3436cfSJed Brown {
33231c3436cfSJed Brown   PetscErrorCode ierr;
33241c3436cfSJed Brown 
33251c3436cfSJed Brown   PetscFunctionBegin;
3326c5033834SJed Brown   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
33271c3436cfSJed Brown   if (vatol) {
33281c3436cfSJed Brown     ierr = PetscObjectReference((PetscObject)vatol);CHKERRQ(ierr);
33291c3436cfSJed Brown     ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
33301c3436cfSJed Brown     ts->vatol = vatol;
33311c3436cfSJed Brown   }
3332c5033834SJed Brown   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
33331c3436cfSJed Brown   if (vrtol) {
33341c3436cfSJed Brown     ierr = PetscObjectReference((PetscObject)vrtol);CHKERRQ(ierr);
33351c3436cfSJed Brown     ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
33361c3436cfSJed Brown     ts->vrtol = vrtol;
33371c3436cfSJed Brown   }
33381c3436cfSJed Brown   PetscFunctionReturn(0);
33391c3436cfSJed Brown }
33401c3436cfSJed Brown 
33411c3436cfSJed Brown #undef __FUNCT__
3342c5033834SJed Brown #define __FUNCT__ "TSGetTolerances"
3343c5033834SJed Brown /*@
3344c5033834SJed Brown    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
3345c5033834SJed Brown 
3346c5033834SJed Brown    Logically Collective
3347c5033834SJed Brown 
3348c5033834SJed Brown    Input Arguments:
3349c5033834SJed Brown .  ts - time integration context
3350c5033834SJed Brown 
3351c5033834SJed Brown    Output Arguments:
3352c5033834SJed Brown +  atol - scalar absolute tolerances, PETSC_NULL to ignore
3353c5033834SJed Brown .  vatol - vector of absolute tolerances, PETSC_NULL to ignore
3354c5033834SJed Brown .  rtol - scalar relative tolerances, PETSC_NULL to ignore
3355c5033834SJed Brown -  vrtol - vector of relative tolerances, PETSC_NULL to ignore
3356c5033834SJed Brown 
3357c5033834SJed Brown    Level: beginner
3358c5033834SJed Brown 
3359c5033834SJed Brown .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
3360c5033834SJed Brown @*/
3361c5033834SJed Brown PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
3362c5033834SJed Brown {
3363c5033834SJed Brown 
3364c5033834SJed Brown   PetscFunctionBegin;
3365c5033834SJed Brown   if (atol)  *atol  = ts->atol;
3366c5033834SJed Brown   if (vatol) *vatol = ts->vatol;
3367c5033834SJed Brown   if (rtol)  *rtol  = ts->rtol;
3368c5033834SJed Brown   if (vrtol) *vrtol = ts->vrtol;
3369c5033834SJed Brown   PetscFunctionReturn(0);
3370c5033834SJed Brown }
3371c5033834SJed Brown 
3372c5033834SJed Brown #undef __FUNCT__
33731c3436cfSJed Brown #define __FUNCT__ "TSErrorNormWRMS"
33741c3436cfSJed Brown /*@
33751c3436cfSJed Brown    TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state
33761c3436cfSJed Brown 
33771c3436cfSJed Brown    Collective on TS
33781c3436cfSJed Brown 
33791c3436cfSJed Brown    Input Arguments:
33801c3436cfSJed Brown +  ts - time stepping context
33811c3436cfSJed Brown -  Y - state vector to be compared to ts->vec_sol
33821c3436cfSJed Brown 
33831c3436cfSJed Brown    Output Arguments:
33841c3436cfSJed Brown .  norm - weighted norm, a value of 1.0 is considered small
33851c3436cfSJed Brown 
33861c3436cfSJed Brown    Level: developer
33871c3436cfSJed Brown 
33881c3436cfSJed Brown .seealso: TSSetTolerances()
33891c3436cfSJed Brown @*/
33901c3436cfSJed Brown PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
33911c3436cfSJed Brown {
33921c3436cfSJed Brown   PetscErrorCode ierr;
33931c3436cfSJed Brown   PetscInt i,n,N;
33941c3436cfSJed Brown   const PetscScalar *x,*y;
33951c3436cfSJed Brown   Vec X;
33961c3436cfSJed Brown   PetscReal sum,gsum;
33971c3436cfSJed Brown 
33981c3436cfSJed Brown   PetscFunctionBegin;
33991c3436cfSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
34001c3436cfSJed Brown   PetscValidHeaderSpecific(Y,VEC_CLASSID,2);
34011c3436cfSJed Brown   PetscValidPointer(norm,3);
34021c3436cfSJed Brown   X = ts->vec_sol;
34031c3436cfSJed Brown   PetscCheckSameTypeAndComm(X,1,Y,2);
34041c3436cfSJed Brown   if (X == Y) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector");
34051c3436cfSJed Brown 
34061c3436cfSJed Brown   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
34071c3436cfSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
34081c3436cfSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
34091c3436cfSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
34101c3436cfSJed Brown   sum = 0.;
3411ceefc113SJed Brown   if (ts->vatol && ts->vrtol) {
3412ceefc113SJed Brown     const PetscScalar *atol,*rtol;
3413ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3414ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3415ceefc113SJed Brown     for (i=0; i<n; i++) {
3416ceefc113SJed Brown       PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3417ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3418ceefc113SJed Brown     }
3419ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3420ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3421ceefc113SJed Brown   } else if (ts->vatol) {       /* vector atol, scalar rtol */
3422ceefc113SJed Brown     const PetscScalar *atol;
3423ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3424ceefc113SJed Brown     for (i=0; i<n; i++) {
3425ceefc113SJed Brown       PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3426ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3427ceefc113SJed Brown     }
3428ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3429ceefc113SJed Brown   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
3430ceefc113SJed Brown     const PetscScalar *rtol;
3431ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3432ceefc113SJed Brown     for (i=0; i<n; i++) {
3433ceefc113SJed Brown       PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3434ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3435ceefc113SJed Brown     }
3436ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3437ceefc113SJed Brown   } else {                      /* scalar atol, scalar rtol */
34381c3436cfSJed Brown     for (i=0; i<n; i++) {
34391c3436cfSJed Brown       PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
34401c3436cfSJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
34411c3436cfSJed Brown     }
3442ceefc113SJed Brown   }
34431c3436cfSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
34441c3436cfSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
34451c3436cfSJed Brown 
34461c3436cfSJed Brown   ierr = MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);CHKERRQ(ierr);
34471c3436cfSJed Brown   *norm = PetscSqrtReal(gsum / N);
34481c3436cfSJed Brown   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
34491c3436cfSJed Brown   PetscFunctionReturn(0);
34501c3436cfSJed Brown }
34511c3436cfSJed Brown 
34521c3436cfSJed Brown #undef __FUNCT__
34538d59e960SJed Brown #define __FUNCT__ "TSSetCFLTimeLocal"
34548d59e960SJed Brown /*@
34558d59e960SJed Brown    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
34568d59e960SJed Brown 
34578d59e960SJed Brown    Logically Collective on TS
34588d59e960SJed Brown 
34598d59e960SJed Brown    Input Arguments:
34608d59e960SJed Brown +  ts - time stepping context
34618d59e960SJed Brown -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)
34628d59e960SJed Brown 
34638d59e960SJed Brown    Note:
34648d59e960SJed Brown    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
34658d59e960SJed Brown 
34668d59e960SJed Brown    Level: intermediate
34678d59e960SJed Brown 
34688d59e960SJed Brown .seealso: TSGetCFLTime(), TSADAPTCFL
34698d59e960SJed Brown @*/
34708d59e960SJed Brown PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
34718d59e960SJed Brown {
34728d59e960SJed Brown 
34738d59e960SJed Brown   PetscFunctionBegin;
34748d59e960SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
34758d59e960SJed Brown   ts->cfltime_local = cfltime;
34768d59e960SJed Brown   ts->cfltime = -1.;
34778d59e960SJed Brown   PetscFunctionReturn(0);
34788d59e960SJed Brown }
34798d59e960SJed Brown 
34808d59e960SJed Brown #undef __FUNCT__
34818d59e960SJed Brown #define __FUNCT__ "TSGetCFLTime"
34828d59e960SJed Brown /*@
34838d59e960SJed Brown    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
34848d59e960SJed Brown 
34858d59e960SJed Brown    Collective on TS
34868d59e960SJed Brown 
34878d59e960SJed Brown    Input Arguments:
34888d59e960SJed Brown .  ts - time stepping context
34898d59e960SJed Brown 
34908d59e960SJed Brown    Output Arguments:
34918d59e960SJed Brown .  cfltime - maximum stable time step for forward Euler
34928d59e960SJed Brown 
34938d59e960SJed Brown    Level: advanced
34948d59e960SJed Brown 
34958d59e960SJed Brown .seealso: TSSetCFLTimeLocal()
34968d59e960SJed Brown @*/
34978d59e960SJed Brown PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
34988d59e960SJed Brown {
34998d59e960SJed Brown   PetscErrorCode ierr;
35008d59e960SJed Brown 
35018d59e960SJed Brown   PetscFunctionBegin;
35028d59e960SJed Brown   if (ts->cfltime < 0) {
35038d59e960SJed Brown     ierr = MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr);
35048d59e960SJed Brown   }
35058d59e960SJed Brown   *cfltime = ts->cfltime;
35068d59e960SJed Brown   PetscFunctionReturn(0);
35078d59e960SJed Brown }
35088d59e960SJed Brown 
35098d59e960SJed Brown #undef __FUNCT__
3510d6ebe24aSShri Abhyankar #define __FUNCT__ "TSVISetVariableBounds"
3511d6ebe24aSShri Abhyankar /*@
3512d6ebe24aSShri Abhyankar    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
3513d6ebe24aSShri Abhyankar 
3514d6ebe24aSShri Abhyankar    Input Parameters:
3515d6ebe24aSShri Abhyankar .  ts   - the TS context.
3516d6ebe24aSShri Abhyankar .  xl   - lower bound.
3517d6ebe24aSShri Abhyankar .  xu   - upper bound.
3518d6ebe24aSShri Abhyankar 
3519d6ebe24aSShri Abhyankar    Notes:
3520d6ebe24aSShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
352133c0c2ddSJed Brown    SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp().
3522d6ebe24aSShri Abhyankar 
35232bd2b0e6SSatish Balay    Level: advanced
35242bd2b0e6SSatish Balay 
3525d6ebe24aSShri Abhyankar @*/
3526d6ebe24aSShri Abhyankar PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
3527d6ebe24aSShri Abhyankar {
3528d6ebe24aSShri Abhyankar   PetscErrorCode ierr;
3529d6ebe24aSShri Abhyankar   SNES           snes;
3530d6ebe24aSShri Abhyankar 
3531d6ebe24aSShri Abhyankar   PetscFunctionBegin;
3532d6ebe24aSShri Abhyankar   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3533d6ebe24aSShri Abhyankar   ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr);
3534d6ebe24aSShri Abhyankar   PetscFunctionReturn(0);
3535d6ebe24aSShri Abhyankar }
3536d6ebe24aSShri Abhyankar 
3537325fc9f4SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
3538c6db04a5SJed Brown #include <mex.h>
3539325fc9f4SBarry Smith 
3540325fc9f4SBarry Smith typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
3541325fc9f4SBarry Smith 
3542325fc9f4SBarry Smith #undef __FUNCT__
3543325fc9f4SBarry Smith #define __FUNCT__ "TSComputeFunction_Matlab"
3544325fc9f4SBarry Smith /*
3545325fc9f4SBarry Smith    TSComputeFunction_Matlab - Calls the function that has been set with
3546325fc9f4SBarry Smith                          TSSetFunctionMatlab().
3547325fc9f4SBarry Smith 
3548325fc9f4SBarry Smith    Collective on TS
3549325fc9f4SBarry Smith 
3550325fc9f4SBarry Smith    Input Parameters:
3551325fc9f4SBarry Smith +  snes - the TS context
3552325fc9f4SBarry Smith -  x - input vector
3553325fc9f4SBarry Smith 
3554325fc9f4SBarry Smith    Output Parameter:
3555325fc9f4SBarry Smith .  y - function vector, as set by TSSetFunction()
3556325fc9f4SBarry Smith 
3557325fc9f4SBarry Smith    Notes:
3558325fc9f4SBarry Smith    TSComputeFunction() is typically used within nonlinear solvers
3559325fc9f4SBarry Smith    implementations, so most users would not generally call this routine
3560325fc9f4SBarry Smith    themselves.
3561325fc9f4SBarry Smith 
3562325fc9f4SBarry Smith    Level: developer
3563325fc9f4SBarry Smith 
3564325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
3565325fc9f4SBarry Smith 
3566325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3567325fc9f4SBarry Smith */
35687087cfbeSBarry Smith PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
3569325fc9f4SBarry Smith {
3570325fc9f4SBarry Smith   PetscErrorCode   ierr;
3571325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3572325fc9f4SBarry Smith   int              nlhs = 1,nrhs = 7;
3573325fc9f4SBarry Smith   mxArray          *plhs[1],*prhs[7];
3574325fc9f4SBarry Smith   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
3575325fc9f4SBarry Smith 
3576325fc9f4SBarry Smith   PetscFunctionBegin;
3577325fc9f4SBarry Smith   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
3578325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3579325fc9f4SBarry Smith   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
3580325fc9f4SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
3581325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,x,3);
3582325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,y,5);
3583325fc9f4SBarry Smith 
3584325fc9f4SBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3585325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3586880f3077SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
3587325fc9f4SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
3588325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3589325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar(time);
3590325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
3591325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3592325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)ly);
3593325fc9f4SBarry Smith   prhs[5] =  mxCreateString(sctx->funcname);
3594325fc9f4SBarry Smith   prhs[6] =  sctx->ctx;
3595325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
3596325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3597325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
3598325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
3599325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
3600325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
3601325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
3602325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
3603325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
3604325fc9f4SBarry Smith   PetscFunctionReturn(0);
3605325fc9f4SBarry Smith }
3606325fc9f4SBarry Smith 
3607325fc9f4SBarry Smith 
3608325fc9f4SBarry Smith #undef __FUNCT__
3609325fc9f4SBarry Smith #define __FUNCT__ "TSSetFunctionMatlab"
3610325fc9f4SBarry Smith /*
3611325fc9f4SBarry Smith    TSSetFunctionMatlab - Sets the function evaluation routine and function
3612325fc9f4SBarry Smith    vector for use by the TS routines in solving ODEs
3613e3c5b3baSBarry Smith    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
3614325fc9f4SBarry Smith 
3615325fc9f4SBarry Smith    Logically Collective on TS
3616325fc9f4SBarry Smith 
3617325fc9f4SBarry Smith    Input Parameters:
3618325fc9f4SBarry Smith +  ts - the TS context
3619325fc9f4SBarry Smith -  func - function evaluation routine
3620325fc9f4SBarry Smith 
3621325fc9f4SBarry Smith    Calling sequence of func:
3622325fc9f4SBarry Smith $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
3623325fc9f4SBarry Smith 
3624325fc9f4SBarry Smith    Level: beginner
3625325fc9f4SBarry Smith 
3626325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
3627325fc9f4SBarry Smith 
3628325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3629325fc9f4SBarry Smith */
3630cdcf91faSSean Farley PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
3631325fc9f4SBarry Smith {
3632325fc9f4SBarry Smith   PetscErrorCode  ierr;
3633325fc9f4SBarry Smith   TSMatlabContext *sctx;
3634325fc9f4SBarry Smith 
3635325fc9f4SBarry Smith   PetscFunctionBegin;
3636325fc9f4SBarry Smith   /* currently sctx is memory bleed */
3637325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3638325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3639325fc9f4SBarry Smith   /*
3640325fc9f4SBarry Smith      This should work, but it doesn't
3641325fc9f4SBarry Smith   sctx->ctx = ctx;
3642325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
3643325fc9f4SBarry Smith   */
3644325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
3645cdcf91faSSean Farley   ierr = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
3646325fc9f4SBarry Smith   PetscFunctionReturn(0);
3647325fc9f4SBarry Smith }
3648325fc9f4SBarry Smith 
3649325fc9f4SBarry Smith #undef __FUNCT__
3650325fc9f4SBarry Smith #define __FUNCT__ "TSComputeJacobian_Matlab"
3651325fc9f4SBarry Smith /*
3652325fc9f4SBarry Smith    TSComputeJacobian_Matlab - Calls the function that has been set with
3653325fc9f4SBarry Smith                          TSSetJacobianMatlab().
3654325fc9f4SBarry Smith 
3655325fc9f4SBarry Smith    Collective on TS
3656325fc9f4SBarry Smith 
3657325fc9f4SBarry Smith    Input Parameters:
3658cdcf91faSSean Farley +  ts - the TS context
3659325fc9f4SBarry Smith .  x - input vector
3660325fc9f4SBarry Smith .  A, B - the matrices
3661325fc9f4SBarry Smith -  ctx - user context
3662325fc9f4SBarry Smith 
3663325fc9f4SBarry Smith    Output Parameter:
3664325fc9f4SBarry Smith .  flag - structure of the matrix
3665325fc9f4SBarry Smith 
3666325fc9f4SBarry Smith    Level: developer
3667325fc9f4SBarry Smith 
3668325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
3669325fc9f4SBarry Smith 
3670325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3671325fc9f4SBarry Smith @*/
3672cdcf91faSSean Farley PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
3673325fc9f4SBarry Smith {
3674325fc9f4SBarry Smith   PetscErrorCode  ierr;
3675325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3676325fc9f4SBarry Smith   int             nlhs = 2,nrhs = 9;
3677325fc9f4SBarry Smith   mxArray         *plhs[2],*prhs[9];
3678325fc9f4SBarry Smith   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
3679325fc9f4SBarry Smith 
3680325fc9f4SBarry Smith   PetscFunctionBegin;
3681cdcf91faSSean Farley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3682325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3683325fc9f4SBarry Smith 
3684325fc9f4SBarry Smith   /* call Matlab function in ctx with arguments x and y */
3685325fc9f4SBarry Smith 
3686cdcf91faSSean Farley   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
3687325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3688325fc9f4SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
3689325fc9f4SBarry Smith   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
3690325fc9f4SBarry Smith   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
3691325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3692325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)time);
3693325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
3694325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3695325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)shift);
3696325fc9f4SBarry Smith   prhs[5] =  mxCreateDoubleScalar((double)lA);
3697325fc9f4SBarry Smith   prhs[6] =  mxCreateDoubleScalar((double)lB);
3698325fc9f4SBarry Smith   prhs[7] =  mxCreateString(sctx->funcname);
3699325fc9f4SBarry Smith   prhs[8] =  sctx->ctx;
3700325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
3701325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3702325fc9f4SBarry Smith   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
3703325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
3704325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
3705325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
3706325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
3707325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
3708325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
3709325fc9f4SBarry Smith   mxDestroyArray(prhs[6]);
3710325fc9f4SBarry Smith   mxDestroyArray(prhs[7]);
3711325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
3712325fc9f4SBarry Smith   mxDestroyArray(plhs[1]);
3713325fc9f4SBarry Smith   PetscFunctionReturn(0);
3714325fc9f4SBarry Smith }
3715325fc9f4SBarry Smith 
3716325fc9f4SBarry Smith 
3717325fc9f4SBarry Smith #undef __FUNCT__
3718325fc9f4SBarry Smith #define __FUNCT__ "TSSetJacobianMatlab"
3719325fc9f4SBarry Smith /*
3720325fc9f4SBarry Smith    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
3721e3c5b3baSBarry 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
3722325fc9f4SBarry Smith 
3723325fc9f4SBarry Smith    Logically Collective on TS
3724325fc9f4SBarry Smith 
3725325fc9f4SBarry Smith    Input Parameters:
3726cdcf91faSSean Farley +  ts - the TS context
3727325fc9f4SBarry Smith .  A,B - Jacobian matrices
3728325fc9f4SBarry Smith .  func - function evaluation routine
3729325fc9f4SBarry Smith -  ctx - user context
3730325fc9f4SBarry Smith 
3731325fc9f4SBarry Smith    Calling sequence of func:
3732cdcf91faSSean Farley $    flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
3733325fc9f4SBarry Smith 
3734325fc9f4SBarry Smith 
3735325fc9f4SBarry Smith    Level: developer
3736325fc9f4SBarry Smith 
3737325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
3738325fc9f4SBarry Smith 
3739325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3740325fc9f4SBarry Smith */
3741cdcf91faSSean Farley PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
3742325fc9f4SBarry Smith {
3743325fc9f4SBarry Smith   PetscErrorCode    ierr;
3744325fc9f4SBarry Smith   TSMatlabContext *sctx;
3745325fc9f4SBarry Smith 
3746325fc9f4SBarry Smith   PetscFunctionBegin;
3747325fc9f4SBarry Smith   /* currently sctx is memory bleed */
3748325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3749325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3750325fc9f4SBarry Smith   /*
3751325fc9f4SBarry Smith      This should work, but it doesn't
3752325fc9f4SBarry Smith   sctx->ctx = ctx;
3753325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
3754325fc9f4SBarry Smith   */
3755325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
3756cdcf91faSSean Farley   ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
3757325fc9f4SBarry Smith   PetscFunctionReturn(0);
3758325fc9f4SBarry Smith }
3759325fc9f4SBarry Smith 
3760b5b1a830SBarry Smith #undef __FUNCT__
3761b5b1a830SBarry Smith #define __FUNCT__ "TSMonitor_Matlab"
3762b5b1a830SBarry Smith /*
3763b5b1a830SBarry Smith    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
3764b5b1a830SBarry Smith 
3765b5b1a830SBarry Smith    Collective on TS
3766b5b1a830SBarry Smith 
3767b5b1a830SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3768b5b1a830SBarry Smith @*/
3769cdcf91faSSean Farley PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
3770b5b1a830SBarry Smith {
3771b5b1a830SBarry Smith   PetscErrorCode  ierr;
3772b5b1a830SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3773a530c242SBarry Smith   int             nlhs = 1,nrhs = 6;
3774b5b1a830SBarry Smith   mxArray         *plhs[1],*prhs[6];
3775b5b1a830SBarry Smith   long long int   lx = 0,ls = 0;
3776b5b1a830SBarry Smith 
3777b5b1a830SBarry Smith   PetscFunctionBegin;
3778cdcf91faSSean Farley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3779b5b1a830SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
3780b5b1a830SBarry Smith 
3781cdcf91faSSean Farley   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
3782b5b1a830SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3783b5b1a830SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3784b5b1a830SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)it);
3785b5b1a830SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)time);
3786b5b1a830SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lx);
3787b5b1a830SBarry Smith   prhs[4] =  mxCreateString(sctx->funcname);
3788b5b1a830SBarry Smith   prhs[5] =  sctx->ctx;
3789b5b1a830SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
3790b5b1a830SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3791b5b1a830SBarry Smith   mxDestroyArray(prhs[0]);
3792b5b1a830SBarry Smith   mxDestroyArray(prhs[1]);
3793b5b1a830SBarry Smith   mxDestroyArray(prhs[2]);
3794b5b1a830SBarry Smith   mxDestroyArray(prhs[3]);
3795b5b1a830SBarry Smith   mxDestroyArray(prhs[4]);
3796b5b1a830SBarry Smith   mxDestroyArray(plhs[0]);
3797b5b1a830SBarry Smith   PetscFunctionReturn(0);
3798b5b1a830SBarry Smith }
3799b5b1a830SBarry Smith 
3800b5b1a830SBarry Smith 
3801b5b1a830SBarry Smith #undef __FUNCT__
3802b5b1a830SBarry Smith #define __FUNCT__ "TSMonitorSetMatlab"
3803b5b1a830SBarry Smith /*
3804b5b1a830SBarry Smith    TSMonitorSetMatlab - Sets the monitor function from Matlab
3805b5b1a830SBarry Smith 
3806b5b1a830SBarry Smith    Level: developer
3807b5b1a830SBarry Smith 
3808b5b1a830SBarry Smith .keywords: TS, nonlinear, set, function
3809b5b1a830SBarry Smith 
3810b5b1a830SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3811b5b1a830SBarry Smith */
3812cdcf91faSSean Farley PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
3813b5b1a830SBarry Smith {
3814b5b1a830SBarry Smith   PetscErrorCode    ierr;
3815b5b1a830SBarry Smith   TSMatlabContext *sctx;
3816b5b1a830SBarry Smith 
3817b5b1a830SBarry Smith   PetscFunctionBegin;
3818b5b1a830SBarry Smith   /* currently sctx is memory bleed */
3819b5b1a830SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3820b5b1a830SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3821b5b1a830SBarry Smith   /*
3822b5b1a830SBarry Smith      This should work, but it doesn't
3823b5b1a830SBarry Smith   sctx->ctx = ctx;
3824b5b1a830SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
3825b5b1a830SBarry Smith   */
3826b5b1a830SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
3827cdcf91faSSean Farley   ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
3828b5b1a830SBarry Smith   PetscFunctionReturn(0);
3829b5b1a830SBarry Smith }
3830325fc9f4SBarry Smith #endif
3831