xref: /petsc/src/ts/interface/ts.c (revision ef20d0600dad68b9ae6c393b3fe871e6fa6b08be)
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) {
117a80ad3e0SBarry Smith       PetscDrawLG lg;
118a80ad3e0SBarry Smith 
119a80ad3e0SBarry Smith       ierr = TSMonitorLGCreate(((PetscObject)ts)->comm,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
120b3603a34SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorLG,lg,(PetscErrorCode (*)(void**))TSMonitorLGDestroy);CHKERRQ(ierr);
121b3603a34SBarry Smith     }
122b3603a34SBarry Smith     opt  = PETSC_FALSE;
123b3603a34SBarry Smith     ierr = PetscOptionsBool("-ts_monitor_solutionode","Monitor solution graphically","TSMonitorSolutionODE",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
124b3603a34SBarry Smith     if (opt) {
125b3603a34SBarry Smith       PetscDrawLG lg;
126b3603a34SBarry Smith 
127ab352700SBarry Smith       ierr = TSMonitorSolutionODECreate(((PetscObject)ts)->comm,3,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,&lg);
128b3603a34SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolutionODE,lg,(PetscErrorCode (*)(void**))TSMonitorSolutionODEDestroy);CHKERRQ(ierr);
129bdad233fSMatthew Knepley     }
13090d69ab7SBarry Smith     opt  = PETSC_FALSE;
131*ef20d060SBarry Smith     ierr = PetscOptionsBool("-ts_monitor_errorode","Monitor error graphically","TSMonitorErrorODE",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
132*ef20d060SBarry Smith     if (opt) {
133*ef20d060SBarry Smith       PetscDrawLG lg;
134*ef20d060SBarry Smith 
135*ef20d060SBarry Smith       ierr = TSMonitorErrorODECreate(((PetscObject)ts)->comm,3,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,&lg);
136*ef20d060SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorErrorODE,lg,(PetscErrorCode (*)(void**))TSMonitorErrorODEDestroy);CHKERRQ(ierr);
137*ef20d060SBarry Smith     }
138*ef20d060SBarry Smith     opt  = PETSC_FALSE;
139acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
140a7cc72afSBarry Smith     if (opt) {
1416083293cSBarry Smith       void        *ctx;
142a80ad3e0SBarry Smith       PetscViewer viewer;
143a80ad3e0SBarry Smith 
144ab352700SBarry Smith       ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm,0,0,PETSC_DECIDE,PETSC_DECIDE,600,400,&viewer);
145a80ad3e0SBarry Smith       ierr = TSMonitorSolutionCreate(ts,viewer,&ctx);CHKERRQ(ierr);
146a80ad3e0SBarry Smith       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1476083293cSBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);CHKERRQ(ierr);
148bdad233fSMatthew Knepley     }
149fb1732b5SBarry Smith     opt  = PETSC_FALSE;
150fb1732b5SBarry Smith     ierr = PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
151fb1732b5SBarry Smith     if (flg) {
152fb1732b5SBarry Smith       PetscViewer ctx;
153fb1732b5SBarry Smith       if (monfilename[0]) {
154fb1732b5SBarry Smith         ierr = PetscViewerBinaryOpen(((PetscObject)ts)->comm,monfilename,FILE_MODE_WRITE,&ctx);CHKERRQ(ierr);
155fb1732b5SBarry Smith       } else {
156fb1732b5SBarry Smith         ctx = PETSC_VIEWER_BINARY_(((PetscObject)ts)->comm);
157fb1732b5SBarry Smith       }
158fb1732b5SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
159fb1732b5SBarry Smith     }
160ed81e22dSJed Brown     opt  = PETSC_FALSE;
161ed81e22dSJed 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);
162ed81e22dSJed Brown     if (flg) {
163ed81e22dSJed Brown       const char *ptr,*ptr2;
164ed81e22dSJed Brown       char *filetemplate;
165ed81e22dSJed Brown       if (!monfilename[0]) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
166ed81e22dSJed Brown       /* Do some cursory validation of the input. */
167ed81e22dSJed Brown       ierr = PetscStrstr(monfilename,"%",(char**)&ptr);CHKERRQ(ierr);
168ed81e22dSJed Brown       if (!ptr) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
169ed81e22dSJed Brown       for (ptr++ ; ptr && *ptr; ptr++) {
170ed81e22dSJed Brown         ierr = PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);CHKERRQ(ierr);
171ed81e22dSJed 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");
172ed81e22dSJed Brown         if (ptr2) break;
173ed81e22dSJed Brown       }
174ed81e22dSJed Brown       ierr = PetscStrallocpy(monfilename,&filetemplate);CHKERRQ(ierr);
175ed81e22dSJed Brown       ierr = TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);CHKERRQ(ierr);
176ed81e22dSJed Brown     }
177bdad233fSMatthew Knepley 
1781c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
1791c3436cfSJed Brown     ierr = TSAdaptSetFromOptions(adapt);CHKERRQ(ierr);
1801c3436cfSJed Brown 
181d52bd9f3SBarry Smith     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
182d52bd9f3SBarry Smith     if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
183d52bd9f3SBarry Smith 
184bdad233fSMatthew Knepley     /* Handle specific TS options */
185abc0a331SBarry Smith     if (ts->ops->setfromoptions) {
186bdad233fSMatthew Knepley       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
187bdad233fSMatthew Knepley     }
1885d973c19SBarry Smith 
1895d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
1905d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
191bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
192bdad233fSMatthew Knepley   PetscFunctionReturn(0);
193bdad233fSMatthew Knepley }
194bdad233fSMatthew Knepley 
195bdad233fSMatthew Knepley #undef __FUNCT__
196cdcf91faSSean Farley #undef __FUNCT__
1974a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
198a7a1495cSBarry Smith /*@
1998c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
200a7a1495cSBarry Smith       set with TSSetRHSJacobian().
201a7a1495cSBarry Smith 
202a7a1495cSBarry Smith    Collective on TS and Vec
203a7a1495cSBarry Smith 
204a7a1495cSBarry Smith    Input Parameters:
205316643e7SJed Brown +  ts - the TS context
206a7a1495cSBarry Smith .  t - current timestep
207a7a1495cSBarry Smith -  x - input vector
208a7a1495cSBarry Smith 
209a7a1495cSBarry Smith    Output Parameters:
210a7a1495cSBarry Smith +  A - Jacobian matrix
211a7a1495cSBarry Smith .  B - optional preconditioning matrix
212a7a1495cSBarry Smith -  flag - flag indicating matrix structure
213a7a1495cSBarry Smith 
214a7a1495cSBarry Smith    Notes:
215a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
216a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
217a7a1495cSBarry Smith 
21894b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
219a7a1495cSBarry Smith    flag parameter.
220a7a1495cSBarry Smith 
221a7a1495cSBarry Smith    Level: developer
222a7a1495cSBarry Smith 
223a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
224a7a1495cSBarry Smith 
22594b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
226a7a1495cSBarry Smith @*/
2277087cfbeSBarry Smith PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
228a7a1495cSBarry Smith {
229dfbe8321SBarry Smith   PetscErrorCode ierr;
2300e4ef248SJed Brown   PetscInt       Xstate;
23124989b8cSPeter Brune   DM             dm;
23224989b8cSPeter Brune   TSDM           tsdm;
23324989b8cSPeter Brune   TSRHSJacobian  rhsjacobianfunc;
23424989b8cSPeter Brune   void           *ctx;
23524989b8cSPeter Brune   TSIJacobian    ijacobianfunc;
236a7a1495cSBarry Smith 
237a7a1495cSBarry Smith   PetscFunctionBegin;
2380700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2390700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
240c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
24124989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
24224989b8cSPeter Brune   ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr);
24324989b8cSPeter Brune   ierr = DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);CHKERRQ(ierr);
24424989b8cSPeter Brune   ierr = DMTSGetIJacobian(dm,&ijacobianfunc,PETSC_NULL);CHKERRQ(ierr);
2450e4ef248SJed Brown   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
2460e4ef248SJed Brown   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) {
2470e4ef248SJed Brown     *flg = ts->rhsjacobian.mstructure;
2480e4ef248SJed Brown     PetscFunctionReturn(0);
249f8ede8e7SJed Brown   }
250d90be118SSean Farley 
25124989b8cSPeter Brune   if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
252d90be118SSean Farley 
25324989b8cSPeter Brune   if (rhsjacobianfunc) {
254d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
255a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
256a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
25724989b8cSPeter Brune     ierr = (*rhsjacobianfunc)(ts,t,X,A,B,flg,ctx);CHKERRQ(ierr);
258a7a1495cSBarry Smith     PetscStackPop;
259d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
260a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2610700a824SBarry Smith     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
2620700a824SBarry Smith     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
263ef66eb69SBarry Smith   } else {
264214bc6a2SJed Brown     ierr = MatZeroEntries(*A);CHKERRQ(ierr);
265214bc6a2SJed Brown     if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);}
266214bc6a2SJed Brown     *flg = SAME_NONZERO_PATTERN;
267ef66eb69SBarry Smith   }
2680e4ef248SJed Brown   ts->rhsjacobian.time = t;
2690e4ef248SJed Brown   ts->rhsjacobian.X = X;
2700e4ef248SJed Brown   ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr);
2710e4ef248SJed Brown   ts->rhsjacobian.mstructure = *flg;
272a7a1495cSBarry Smith   PetscFunctionReturn(0);
273a7a1495cSBarry Smith }
274a7a1495cSBarry Smith 
2754a2ae208SSatish Balay #undef __FUNCT__
2764a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
277316643e7SJed Brown /*@
278d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
279d763cef2SBarry Smith 
280316643e7SJed Brown    Collective on TS and Vec
281316643e7SJed Brown 
282316643e7SJed Brown    Input Parameters:
283316643e7SJed Brown +  ts - the TS context
284316643e7SJed Brown .  t - current time
285316643e7SJed Brown -  x - state vector
286316643e7SJed Brown 
287316643e7SJed Brown    Output Parameter:
288316643e7SJed Brown .  y - right hand side
289316643e7SJed Brown 
290316643e7SJed Brown    Note:
291316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
292316643e7SJed Brown    is used internally within the nonlinear solvers.
293316643e7SJed Brown 
294316643e7SJed Brown    Level: developer
295316643e7SJed Brown 
296316643e7SJed Brown .keywords: TS, compute
297316643e7SJed Brown 
298316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction()
299316643e7SJed Brown @*/
300dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
301d763cef2SBarry Smith {
302dfbe8321SBarry Smith   PetscErrorCode ierr;
303d763cef2SBarry Smith 
30424989b8cSPeter Brune   TSRHSFunction rhsfunction;
30524989b8cSPeter Brune   TSIFunction   ifunction;
30624989b8cSPeter Brune   void          *ctx;
30724989b8cSPeter Brune   DM            dm;
30824989b8cSPeter Brune 
309d763cef2SBarry Smith   PetscFunctionBegin;
3100700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3110700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3120700a824SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
31324989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
31424989b8cSPeter Brune   ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctx);CHKERRQ(ierr);
31524989b8cSPeter Brune   ierr = DMTSGetIFunction(dm,&ifunction,PETSC_NULL);CHKERRQ(ierr);
316d763cef2SBarry Smith 
31724989b8cSPeter Brune   if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
318d763cef2SBarry Smith 
319d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
32024989b8cSPeter Brune   if (rhsfunction) {
321d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
32224989b8cSPeter Brune     ierr = (*rhsfunction)(ts,t,x,y,ctx);CHKERRQ(ierr);
323d763cef2SBarry Smith     PetscStackPop;
324214bc6a2SJed Brown   } else {
325214bc6a2SJed Brown     ierr = VecZeroEntries(y);CHKERRQ(ierr);
326b2cd27e8SJed Brown   }
32744a41b28SSean Farley 
328d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
329d763cef2SBarry Smith   PetscFunctionReturn(0);
330d763cef2SBarry Smith }
331d763cef2SBarry Smith 
3324a2ae208SSatish Balay #undef __FUNCT__
333*ef20d060SBarry Smith #define __FUNCT__ "TSComputeSolutionFunction"
334*ef20d060SBarry Smith /*@
335*ef20d060SBarry Smith    TSComputeSolutionFunction - Evaluates the solution function.
336*ef20d060SBarry Smith 
337*ef20d060SBarry Smith    Collective on TS and Vec
338*ef20d060SBarry Smith 
339*ef20d060SBarry Smith    Input Parameters:
340*ef20d060SBarry Smith +  ts - the TS context
341*ef20d060SBarry Smith -  t - current time
342*ef20d060SBarry Smith 
343*ef20d060SBarry Smith    Output Parameter:
344*ef20d060SBarry Smith .  y - right hand side
345*ef20d060SBarry Smith 
346*ef20d060SBarry Smith    Note:
347*ef20d060SBarry Smith    Most users should not need to explicitly call this routine, as it
348*ef20d060SBarry Smith    is used internally within the nonlinear solvers.
349*ef20d060SBarry Smith 
350*ef20d060SBarry Smith    Level: developer
351*ef20d060SBarry Smith 
352*ef20d060SBarry Smith .keywords: TS, compute
353*ef20d060SBarry Smith 
354*ef20d060SBarry Smith .seealso: TSSetRHSFunction(), TSComputeIFunction()
355*ef20d060SBarry Smith @*/
356*ef20d060SBarry Smith PetscErrorCode TSComputeSolutionFunction(TS ts,PetscReal t,Vec x)
357*ef20d060SBarry Smith {
358*ef20d060SBarry Smith   PetscErrorCode     ierr;
359*ef20d060SBarry Smith   TSSolutionFunction solutionfunction;
360*ef20d060SBarry Smith   void               *ctx;
361*ef20d060SBarry Smith   DM                 dm;
362*ef20d060SBarry Smith 
363*ef20d060SBarry Smith   PetscFunctionBegin;
364*ef20d060SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
365*ef20d060SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
366*ef20d060SBarry Smith   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
367*ef20d060SBarry Smith   ierr = DMTSGetSolutionFunction(dm,&solutionfunction,&ctx);CHKERRQ(ierr);
368*ef20d060SBarry Smith 
369*ef20d060SBarry Smith   if (solutionfunction) {
370*ef20d060SBarry Smith     PetscStackPush("TS user right-hand-side function");
371*ef20d060SBarry Smith     ierr = (*solutionfunction)(ts,t,x,ctx);CHKERRQ(ierr);
372*ef20d060SBarry Smith     PetscStackPop;
373*ef20d060SBarry Smith   }
374*ef20d060SBarry Smith   PetscFunctionReturn(0);
375*ef20d060SBarry Smith }
376*ef20d060SBarry Smith 
377*ef20d060SBarry Smith #undef __FUNCT__
378214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSVec_Private"
379214bc6a2SJed Brown static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
380214bc6a2SJed Brown {
3812dd45cf8SJed Brown   Vec            F;
382214bc6a2SJed Brown   PetscErrorCode ierr;
383214bc6a2SJed Brown 
384214bc6a2SJed Brown   PetscFunctionBegin;
3850da9a4f0SJed Brown   *Frhs = PETSC_NULL;
3862dd45cf8SJed Brown   ierr = TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
387214bc6a2SJed Brown   if (!ts->Frhs) {
3882dd45cf8SJed Brown     ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr);
389214bc6a2SJed Brown   }
390214bc6a2SJed Brown   *Frhs = ts->Frhs;
391214bc6a2SJed Brown   PetscFunctionReturn(0);
392214bc6a2SJed Brown }
393214bc6a2SJed Brown 
394214bc6a2SJed Brown #undef __FUNCT__
395214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSMats_Private"
396214bc6a2SJed Brown static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
397214bc6a2SJed Brown {
398214bc6a2SJed Brown   Mat            A,B;
3992dd45cf8SJed Brown   PetscErrorCode ierr;
400214bc6a2SJed Brown 
401214bc6a2SJed Brown   PetscFunctionBegin;
402214bc6a2SJed Brown   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
403214bc6a2SJed Brown   if (Arhs) {
404214bc6a2SJed Brown     if (!ts->Arhs) {
405214bc6a2SJed Brown       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
406214bc6a2SJed Brown     }
407214bc6a2SJed Brown     *Arhs = ts->Arhs;
408214bc6a2SJed Brown   }
409214bc6a2SJed Brown   if (Brhs) {
410214bc6a2SJed Brown     if (!ts->Brhs) {
411214bc6a2SJed Brown       ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
412214bc6a2SJed Brown     }
413214bc6a2SJed Brown     *Brhs = ts->Brhs;
414214bc6a2SJed Brown   }
415214bc6a2SJed Brown   PetscFunctionReturn(0);
416214bc6a2SJed Brown }
417214bc6a2SJed Brown 
418214bc6a2SJed Brown #undef __FUNCT__
419316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction"
420316643e7SJed Brown /*@
421316643e7SJed Brown    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
422316643e7SJed Brown 
423316643e7SJed Brown    Collective on TS and Vec
424316643e7SJed Brown 
425316643e7SJed Brown    Input Parameters:
426316643e7SJed Brown +  ts - the TS context
427316643e7SJed Brown .  t - current time
428316643e7SJed Brown .  X - state vector
429214bc6a2SJed Brown .  Xdot - time derivative of state vector
430214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
431316643e7SJed Brown 
432316643e7SJed Brown    Output Parameter:
433316643e7SJed Brown .  Y - right hand side
434316643e7SJed Brown 
435316643e7SJed Brown    Note:
436316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
437316643e7SJed Brown    is used internally within the nonlinear solvers.
438316643e7SJed Brown 
439316643e7SJed Brown    If the user did did not write their equations in implicit form, this
440316643e7SJed Brown    function recasts them in implicit form.
441316643e7SJed Brown 
442316643e7SJed Brown    Level: developer
443316643e7SJed Brown 
444316643e7SJed Brown .keywords: TS, compute
445316643e7SJed Brown 
446316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction()
447316643e7SJed Brown @*/
448214bc6a2SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
449316643e7SJed Brown {
450316643e7SJed Brown   PetscErrorCode ierr;
45124989b8cSPeter Brune   TSIFunction    ifunction;
45224989b8cSPeter Brune   TSRHSFunction  rhsfunction;
45324989b8cSPeter Brune   void           *ctx;
45424989b8cSPeter Brune   DM             dm;
455316643e7SJed Brown 
456316643e7SJed Brown   PetscFunctionBegin;
4570700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4580700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
4590700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
4600700a824SBarry Smith   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
461316643e7SJed Brown 
46224989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
46324989b8cSPeter Brune   ierr = DMTSGetIFunction(dm,&ifunction,&ctx);CHKERRQ(ierr);
46424989b8cSPeter Brune   ierr = DMTSGetRHSFunction(dm,&rhsfunction,PETSC_NULL);CHKERRQ(ierr);
46524989b8cSPeter Brune 
46624989b8cSPeter Brune   if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
467d90be118SSean Farley 
468316643e7SJed Brown   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
46924989b8cSPeter Brune   if (ifunction) {
470316643e7SJed Brown     PetscStackPush("TS user implicit function");
47124989b8cSPeter Brune     ierr = (*ifunction)(ts,t,X,Xdot,Y,ctx);CHKERRQ(ierr);
472316643e7SJed Brown     PetscStackPop;
473214bc6a2SJed Brown   }
474214bc6a2SJed Brown   if (imex) {
47524989b8cSPeter Brune     if (!ifunction) {
4762dd45cf8SJed Brown       ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);
4772dd45cf8SJed Brown     }
47824989b8cSPeter Brune   } else if (rhsfunction) {
47924989b8cSPeter Brune     if (ifunction) {
480214bc6a2SJed Brown       Vec Frhs;
481214bc6a2SJed Brown       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
482214bc6a2SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr);
483214bc6a2SJed Brown       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
4842dd45cf8SJed Brown     } else {
4852dd45cf8SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr);
4862dd45cf8SJed Brown       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
487316643e7SJed Brown     }
4884a6899ffSJed Brown   }
489316643e7SJed Brown   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
490316643e7SJed Brown   PetscFunctionReturn(0);
491316643e7SJed Brown }
492316643e7SJed Brown 
493316643e7SJed Brown #undef __FUNCT__
494316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian"
495316643e7SJed Brown /*@
496316643e7SJed Brown    TSComputeIJacobian - Evaluates the Jacobian of the DAE
497316643e7SJed Brown 
498316643e7SJed Brown    Collective on TS and Vec
499316643e7SJed Brown 
500316643e7SJed Brown    Input
501316643e7SJed Brown       Input Parameters:
502316643e7SJed Brown +  ts - the TS context
503316643e7SJed Brown .  t - current timestep
504316643e7SJed Brown .  X - state vector
505316643e7SJed Brown .  Xdot - time derivative of state vector
506214bc6a2SJed Brown .  shift - shift to apply, see note below
507214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
508316643e7SJed Brown 
509316643e7SJed Brown    Output Parameters:
510316643e7SJed Brown +  A - Jacobian matrix
511316643e7SJed Brown .  B - optional preconditioning matrix
512316643e7SJed Brown -  flag - flag indicating matrix structure
513316643e7SJed Brown 
514316643e7SJed Brown    Notes:
515316643e7SJed Brown    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
516316643e7SJed Brown 
517316643e7SJed Brown    dF/dX + shift*dF/dXdot
518316643e7SJed Brown 
519316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
520316643e7SJed Brown    is used internally within the nonlinear solvers.
521316643e7SJed Brown 
522316643e7SJed Brown    Level: developer
523316643e7SJed Brown 
524316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix
525316643e7SJed Brown 
526316643e7SJed Brown .seealso:  TSSetIJacobian()
52763495f91SJed Brown @*/
528214bc6a2SJed Brown PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
529316643e7SJed Brown {
5300026cea9SSean Farley   PetscInt Xstate, Xdotstate;
531316643e7SJed Brown   PetscErrorCode ierr;
53224989b8cSPeter Brune   TSIJacobian    ijacobian;
53324989b8cSPeter Brune   TSRHSJacobian  rhsjacobian;
53424989b8cSPeter Brune   DM             dm;
53524989b8cSPeter Brune   void           *ctx;
536316643e7SJed Brown 
537316643e7SJed Brown   PetscFunctionBegin;
53824989b8cSPeter Brune 
5390700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5400700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
5410700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
542316643e7SJed Brown   PetscValidPointer(A,6);
5430700a824SBarry Smith   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
544316643e7SJed Brown   PetscValidPointer(B,7);
5450700a824SBarry Smith   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
546316643e7SJed Brown   PetscValidPointer(flg,8);
54724989b8cSPeter Brune 
54824989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
54924989b8cSPeter Brune   ierr = DMTSGetIJacobian(dm,&ijacobian,&ctx);CHKERRQ(ierr);
55024989b8cSPeter Brune   ierr = DMTSGetRHSJacobian(dm,&rhsjacobian,PETSC_NULL);CHKERRQ(ierr);
55124989b8cSPeter Brune 
5520026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
5530026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr);
5540026cea9SSean 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))) {
5550026cea9SSean Farley     *flg = ts->ijacobian.mstructure;
5560026cea9SSean Farley     ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr);
5570026cea9SSean Farley     PetscFunctionReturn(0);
5580026cea9SSean Farley   }
559316643e7SJed Brown 
56024989b8cSPeter Brune   if (!rhsjacobian && !ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
561316643e7SJed Brown 
5624e684422SJed Brown   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
563316643e7SJed Brown   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
56424989b8cSPeter Brune   if (ijacobian) {
5652dd45cf8SJed Brown     *flg = DIFFERENT_NONZERO_PATTERN;
566316643e7SJed Brown     PetscStackPush("TS user implicit Jacobian");
56724989b8cSPeter Brune     ierr = (*ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ctx);CHKERRQ(ierr);
568316643e7SJed Brown     PetscStackPop;
569214bc6a2SJed Brown     /* make sure user returned a correct Jacobian and preconditioner */
570214bc6a2SJed Brown     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
571214bc6a2SJed Brown     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
5724a6899ffSJed Brown   }
573214bc6a2SJed Brown   if (imex) {
57424989b8cSPeter Brune     if (!ijacobian) {  /* system was written as Xdot = F(t,X) */
575214bc6a2SJed Brown       ierr = MatZeroEntries(*A);CHKERRQ(ierr);
5762dd45cf8SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
577214bc6a2SJed Brown       if (*A != *B) {
578214bc6a2SJed Brown         ierr = MatZeroEntries(*B);CHKERRQ(ierr);
579214bc6a2SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
580214bc6a2SJed Brown       }
581214bc6a2SJed Brown       *flg = SAME_PRECONDITIONER;
582214bc6a2SJed Brown     }
583214bc6a2SJed Brown   } else {
58424989b8cSPeter Brune     if (!ijacobian) {
585214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr);
586214bc6a2SJed Brown       ierr = MatScale(*A,-1);CHKERRQ(ierr);
587214bc6a2SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
588316643e7SJed Brown       if (*A != *B) {
589316643e7SJed Brown         ierr = MatScale(*B,-1);CHKERRQ(ierr);
590316643e7SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
591316643e7SJed Brown       }
59224989b8cSPeter Brune     } else if (rhsjacobian) {
593214bc6a2SJed Brown       Mat Arhs,Brhs;
594214bc6a2SJed Brown       MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
595214bc6a2SJed Brown       ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
596214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
597214bc6a2SJed Brown       axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
598214bc6a2SJed Brown       ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr);
599214bc6a2SJed Brown       if (*A != *B) {
600214bc6a2SJed Brown         ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr);
601214bc6a2SJed Brown       }
602214bc6a2SJed Brown       *flg = PetscMin(*flg,flg2);
603214bc6a2SJed Brown     }
604316643e7SJed Brown   }
6050026cea9SSean Farley 
6060026cea9SSean Farley   ts->ijacobian.time = t;
6070026cea9SSean Farley   ts->ijacobian.X = X;
6080026cea9SSean Farley   ts->ijacobian.Xdot = Xdot;
6090026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr);
6100026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr);
6110026cea9SSean Farley   ts->ijacobian.shift = shift;
6120026cea9SSean Farley   ts->ijacobian.imex = imex;
6130026cea9SSean Farley   ts->ijacobian.mstructure = *flg;
614316643e7SJed Brown   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
615316643e7SJed Brown   PetscFunctionReturn(0);
616316643e7SJed Brown }
617316643e7SJed Brown 
618316643e7SJed Brown #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
620d763cef2SBarry Smith /*@C
621d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
622d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
623d763cef2SBarry Smith 
6243f9fe445SBarry Smith     Logically Collective on TS
625d763cef2SBarry Smith 
626d763cef2SBarry Smith     Input Parameters:
627d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
628ca94891dSJed Brown .   r - vector to put the computed right hand side (or PETSC_NULL to have it created)
629d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
630d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
631d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
632d763cef2SBarry Smith 
633d763cef2SBarry Smith     Calling sequence of func:
63487828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
635d763cef2SBarry Smith 
636d763cef2SBarry Smith +   t - current timestep
637d763cef2SBarry Smith .   u - input vector
638d763cef2SBarry Smith .   F - function vector
639d763cef2SBarry Smith -   ctx - [optional] user-defined function context
640d763cef2SBarry Smith 
641d763cef2SBarry Smith     Level: beginner
642d763cef2SBarry Smith 
643d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
644d763cef2SBarry Smith 
645d6cbdb99SBarry Smith .seealso: TSSetRHSJacobian(), TSSetIJacobian()
646d763cef2SBarry Smith @*/
647089b2837SJed Brown PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
648d763cef2SBarry Smith {
649089b2837SJed Brown   PetscErrorCode ierr;
650089b2837SJed Brown   SNES           snes;
651e856ceecSJed Brown   Vec            ralloc = PETSC_NULL;
65224989b8cSPeter Brune   DM             dm;
653d763cef2SBarry Smith 
654089b2837SJed Brown   PetscFunctionBegin;
6550700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
656ca94891dSJed Brown   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
65724989b8cSPeter Brune 
65824989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
65924989b8cSPeter Brune   ierr = DMTSSetRHSFunction(dm,f,ctx);CHKERRQ(ierr);
660089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
661e856ceecSJed Brown   if (!r && !ts->dm && ts->vec_sol) {
662e856ceecSJed Brown     ierr = VecDuplicate(ts->vec_sol,&ralloc);CHKERRQ(ierr);
663e856ceecSJed Brown     r = ralloc;
664e856ceecSJed Brown   }
665089b2837SJed Brown   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
666e856ceecSJed Brown   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
667d763cef2SBarry Smith   PetscFunctionReturn(0);
668d763cef2SBarry Smith }
669d763cef2SBarry Smith 
6704a2ae208SSatish Balay #undef __FUNCT__
671*ef20d060SBarry Smith #define __FUNCT__ "TSSetSolutionFunction"
672*ef20d060SBarry Smith /*@C
673*ef20d060SBarry Smith     TSSetSolutionFunction - Provide a function that computes the solution of the ODE
674*ef20d060SBarry Smith 
675*ef20d060SBarry Smith     Logically Collective on TS
676*ef20d060SBarry Smith 
677*ef20d060SBarry Smith     Input Parameters:
678*ef20d060SBarry Smith +   ts - the TS context obtained from TSCreate()
679*ef20d060SBarry Smith .   f - routine for evaluating the solution
680*ef20d060SBarry Smith -   ctx - [optional] user-defined context for private data for the
681*ef20d060SBarry Smith           function evaluation routine (may be PETSC_NULL)
682*ef20d060SBarry Smith 
683*ef20d060SBarry Smith     Calling sequence of func:
684*ef20d060SBarry Smith $     func (TS ts,PetscReal t,Vec u,void *ctx);
685*ef20d060SBarry Smith 
686*ef20d060SBarry Smith +   t - current timestep
687*ef20d060SBarry Smith .   u - output vector
688*ef20d060SBarry Smith -   ctx - [optional] user-defined function context
689*ef20d060SBarry Smith 
690*ef20d060SBarry Smith     Level: beginner
691*ef20d060SBarry Smith 
692*ef20d060SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
693*ef20d060SBarry Smith 
694*ef20d060SBarry Smith .seealso: TSSetRHSJacobian(), TSSetIJacobian()
695*ef20d060SBarry Smith @*/
696*ef20d060SBarry Smith PetscErrorCode  TSSetSolutionFunction(TS ts,PetscErrorCode (*f)(TS,PetscReal,Vec,void*),void *ctx)
697*ef20d060SBarry Smith {
698*ef20d060SBarry Smith   PetscErrorCode ierr;
699*ef20d060SBarry Smith   DM             dm;
700*ef20d060SBarry Smith 
701*ef20d060SBarry Smith   PetscFunctionBegin;
702*ef20d060SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
703*ef20d060SBarry Smith   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
704*ef20d060SBarry Smith   ierr = DMTSSetSolutionFunction(dm,f,ctx);CHKERRQ(ierr);
705*ef20d060SBarry Smith   PetscFunctionReturn(0);
706*ef20d060SBarry Smith }
707*ef20d060SBarry Smith 
708*ef20d060SBarry Smith #undef __FUNCT__
7094a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
710d763cef2SBarry Smith /*@C
711d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
712d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
713d763cef2SBarry Smith 
7143f9fe445SBarry Smith    Logically Collective on TS
715d763cef2SBarry Smith 
716d763cef2SBarry Smith    Input Parameters:
717d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
718d763cef2SBarry Smith .  A   - Jacobian matrix
719d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
720d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
721d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
722d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
723d763cef2SBarry Smith 
724d763cef2SBarry Smith    Calling sequence of func:
72587828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
726d763cef2SBarry Smith 
727d763cef2SBarry Smith +  t - current timestep
728d763cef2SBarry Smith .  u - input vector
729d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
730d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
731d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
73294b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
733d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
734d763cef2SBarry Smith 
735d763cef2SBarry Smith    Notes:
73694b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
737d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
738d763cef2SBarry Smith 
739d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
740d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
74156335db2SHong Zhang    completely new matrix structure (not just different matrix elements)
742d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
743d763cef2SBarry Smith    throughout the global iterations.
744d763cef2SBarry Smith 
745d763cef2SBarry Smith    Level: beginner
746d763cef2SBarry Smith 
747d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
748d763cef2SBarry Smith 
749ab637aeaSJed Brown .seealso: SNESDefaultComputeJacobianColor(), TSSetRHSFunction()
750d763cef2SBarry Smith 
751d763cef2SBarry Smith @*/
752089b2837SJed Brown PetscErrorCode  TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
753d763cef2SBarry Smith {
754277b19d0SLisandro Dalcin   PetscErrorCode ierr;
755089b2837SJed Brown   SNES           snes;
75624989b8cSPeter Brune   DM             dm;
75724989b8cSPeter Brune   TSIJacobian    ijacobian;
758277b19d0SLisandro Dalcin 
759d763cef2SBarry Smith   PetscFunctionBegin;
7600700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
761277b19d0SLisandro Dalcin   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
762277b19d0SLisandro Dalcin   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
763277b19d0SLisandro Dalcin   if (A) PetscCheckSameComm(ts,1,A,2);
764277b19d0SLisandro Dalcin   if (B) PetscCheckSameComm(ts,1,B,3);
765d763cef2SBarry Smith 
76624989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
76724989b8cSPeter Brune   ierr = DMTSSetRHSJacobian(dm,f,ctx);CHKERRQ(ierr);
76824989b8cSPeter Brune   ierr = DMTSGetIJacobian(dm,&ijacobian,PETSC_NULL);CHKERRQ(ierr);
76924989b8cSPeter Brune 
770089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
7715f659677SPeter Brune   if (!ijacobian) {
772089b2837SJed Brown     ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
7730e4ef248SJed Brown   }
7740e4ef248SJed Brown   if (A) {
7750e4ef248SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
7760e4ef248SJed Brown     ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
7770e4ef248SJed Brown     ts->Arhs = A;
7780e4ef248SJed Brown   }
7790e4ef248SJed Brown   if (B) {
7800e4ef248SJed Brown     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
7810e4ef248SJed Brown     ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
7820e4ef248SJed Brown     ts->Brhs = B;
7830e4ef248SJed Brown   }
784d763cef2SBarry Smith   PetscFunctionReturn(0);
785d763cef2SBarry Smith }
786d763cef2SBarry Smith 
787316643e7SJed Brown 
788316643e7SJed Brown #undef __FUNCT__
789316643e7SJed Brown #define __FUNCT__ "TSSetIFunction"
790316643e7SJed Brown /*@C
791316643e7SJed Brown    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
792316643e7SJed Brown 
7933f9fe445SBarry Smith    Logically Collective on TS
794316643e7SJed Brown 
795316643e7SJed Brown    Input Parameters:
796316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
797ca94891dSJed Brown .  r   - vector to hold the residual (or PETSC_NULL to have it created internally)
798316643e7SJed Brown .  f   - the function evaluation routine
799316643e7SJed Brown -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
800316643e7SJed Brown 
801316643e7SJed Brown    Calling sequence of f:
802316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
803316643e7SJed Brown 
804316643e7SJed Brown +  t   - time at step/stage being solved
805316643e7SJed Brown .  u   - state vector
806316643e7SJed Brown .  u_t - time derivative of state vector
807316643e7SJed Brown .  F   - function vector
808316643e7SJed Brown -  ctx - [optional] user-defined context for matrix evaluation routine
809316643e7SJed Brown 
810316643e7SJed Brown    Important:
811d6cbdb99SBarry Smith    The user MUST call either this routine, TSSetRHSFunction().  This routine must be used when not solving an ODE, for example a DAE.
812316643e7SJed Brown 
813316643e7SJed Brown    Level: beginner
814316643e7SJed Brown 
815316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian
816316643e7SJed Brown 
817d6cbdb99SBarry Smith .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
818316643e7SJed Brown @*/
819089b2837SJed Brown PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
820316643e7SJed Brown {
821089b2837SJed Brown   PetscErrorCode ierr;
822089b2837SJed Brown   SNES           snes;
823e856ceecSJed Brown   Vec            resalloc = PETSC_NULL;
82424989b8cSPeter Brune   DM             dm;
825316643e7SJed Brown 
826316643e7SJed Brown   PetscFunctionBegin;
8270700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
828ca94891dSJed Brown   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
82924989b8cSPeter Brune 
83024989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
83124989b8cSPeter Brune   ierr = DMTSSetIFunction(dm,f,ctx);CHKERRQ(ierr);
83224989b8cSPeter Brune 
833089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
834e856ceecSJed Brown   if (!res && !ts->dm && ts->vec_sol) {
835e856ceecSJed Brown     ierr = VecDuplicate(ts->vec_sol,&resalloc);CHKERRQ(ierr);
836e856ceecSJed Brown     res = resalloc;
837e856ceecSJed Brown   }
838089b2837SJed Brown   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
839e856ceecSJed Brown   ierr = VecDestroy(&resalloc);CHKERRQ(ierr);
84024989b8cSPeter Brune 
841089b2837SJed Brown   PetscFunctionReturn(0);
842089b2837SJed Brown }
843089b2837SJed Brown 
844089b2837SJed Brown #undef __FUNCT__
845089b2837SJed Brown #define __FUNCT__ "TSGetIFunction"
846089b2837SJed Brown /*@C
847089b2837SJed Brown    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
848089b2837SJed Brown 
849089b2837SJed Brown    Not Collective
850089b2837SJed Brown 
851089b2837SJed Brown    Input Parameter:
852089b2837SJed Brown .  ts - the TS context
853089b2837SJed Brown 
854089b2837SJed Brown    Output Parameter:
855089b2837SJed Brown +  r - vector to hold residual (or PETSC_NULL)
856089b2837SJed Brown .  func - the function to compute residual (or PETSC_NULL)
857089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
858089b2837SJed Brown 
859089b2837SJed Brown    Level: advanced
860089b2837SJed Brown 
861089b2837SJed Brown .keywords: TS, nonlinear, get, function
862089b2837SJed Brown 
863089b2837SJed Brown .seealso: TSSetIFunction(), SNESGetFunction()
864089b2837SJed Brown @*/
865089b2837SJed Brown PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
866089b2837SJed Brown {
867089b2837SJed Brown   PetscErrorCode ierr;
868089b2837SJed Brown   SNES snes;
86924989b8cSPeter Brune   DM   dm;
870089b2837SJed Brown 
871089b2837SJed Brown   PetscFunctionBegin;
872089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
873089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
874089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
87524989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
87624989b8cSPeter Brune   ierr = DMTSGetIFunction(dm,func,ctx);CHKERRQ(ierr);
877089b2837SJed Brown   PetscFunctionReturn(0);
878089b2837SJed Brown }
879089b2837SJed Brown 
880089b2837SJed Brown #undef __FUNCT__
881089b2837SJed Brown #define __FUNCT__ "TSGetRHSFunction"
882089b2837SJed Brown /*@C
883089b2837SJed Brown    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
884089b2837SJed Brown 
885089b2837SJed Brown    Not Collective
886089b2837SJed Brown 
887089b2837SJed Brown    Input Parameter:
888089b2837SJed Brown .  ts - the TS context
889089b2837SJed Brown 
890089b2837SJed Brown    Output Parameter:
891089b2837SJed Brown +  r - vector to hold computed right hand side (or PETSC_NULL)
892089b2837SJed Brown .  func - the function to compute right hand side (or PETSC_NULL)
893089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
894089b2837SJed Brown 
895089b2837SJed Brown    Level: advanced
896089b2837SJed Brown 
897089b2837SJed Brown .keywords: TS, nonlinear, get, function
898089b2837SJed Brown 
899089b2837SJed Brown .seealso: TSSetRhsfunction(), SNESGetFunction()
900089b2837SJed Brown @*/
901089b2837SJed Brown PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
902089b2837SJed Brown {
903089b2837SJed Brown   PetscErrorCode ierr;
904089b2837SJed Brown   SNES snes;
90524989b8cSPeter Brune   DM   dm;
906089b2837SJed Brown 
907089b2837SJed Brown   PetscFunctionBegin;
908089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
909089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
910089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
91124989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
91224989b8cSPeter Brune   ierr = DMTSGetRHSFunction(dm,func,ctx);CHKERRQ(ierr);
913316643e7SJed Brown   PetscFunctionReturn(0);
914316643e7SJed Brown }
915316643e7SJed Brown 
916316643e7SJed Brown #undef __FUNCT__
917316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian"
918316643e7SJed Brown /*@C
919a4f0a591SBarry Smith    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
920a4f0a591SBarry Smith         you provided with TSSetIFunction().
921316643e7SJed Brown 
9223f9fe445SBarry Smith    Logically Collective on TS
923316643e7SJed Brown 
924316643e7SJed Brown    Input Parameters:
925316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
926316643e7SJed Brown .  A   - Jacobian matrix
927316643e7SJed Brown .  B   - preconditioning matrix for A (may be same as A)
928316643e7SJed Brown .  f   - the Jacobian evaluation routine
929316643e7SJed Brown -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
930316643e7SJed Brown 
931316643e7SJed Brown    Calling sequence of f:
9321b4a444bSJed Brown $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
933316643e7SJed Brown 
934316643e7SJed Brown +  t    - time at step/stage being solved
9351b4a444bSJed Brown .  U    - state vector
9361b4a444bSJed Brown .  U_t  - time derivative of state vector
937316643e7SJed Brown .  a    - shift
9381b4a444bSJed Brown .  A    - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
939316643e7SJed Brown .  B    - preconditioning matrix for A, may be same as A
940316643e7SJed Brown .  flag - flag indicating information about the preconditioner matrix
941316643e7SJed Brown           structure (same as flag in KSPSetOperators())
942316643e7SJed Brown -  ctx  - [optional] user-defined context for matrix evaluation routine
943316643e7SJed Brown 
944316643e7SJed Brown    Notes:
945316643e7SJed Brown    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
946316643e7SJed Brown 
947a4f0a591SBarry Smith    The matrix dF/dU + a*dF/dU_t you provide turns out to be
948a4f0a591SBarry 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.
949a4f0a591SBarry Smith    The time integrator internally approximates U_t by W+a*U where the positive "shift"
950a4f0a591SBarry Smith    a and vector W depend on the integration method, step size, and past states. For example with
951a4f0a591SBarry Smith    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
952a4f0a591SBarry Smith    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
953a4f0a591SBarry Smith 
954316643e7SJed Brown    Level: beginner
955316643e7SJed Brown 
956316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian
957316643e7SJed Brown 
958ab637aeaSJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESDefaultComputeJacobianColor(), SNESDefaultComputeJacobian()
959316643e7SJed Brown 
960316643e7SJed Brown @*/
9617087cfbeSBarry Smith PetscErrorCode  TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
962316643e7SJed Brown {
963316643e7SJed Brown   PetscErrorCode ierr;
964089b2837SJed Brown   SNES           snes;
96524989b8cSPeter Brune   DM             dm;
966316643e7SJed Brown 
967316643e7SJed Brown   PetscFunctionBegin;
9680700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
9690700a824SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
9700700a824SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
971316643e7SJed Brown   if (A) PetscCheckSameComm(ts,1,A,2);
972316643e7SJed Brown   if (B) PetscCheckSameComm(ts,1,B,3);
97324989b8cSPeter Brune 
97424989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
97524989b8cSPeter Brune   ierr = DMTSSetIJacobian(dm,f,ctx);CHKERRQ(ierr);
97624989b8cSPeter Brune 
977089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
978089b2837SJed Brown   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
979316643e7SJed Brown   PetscFunctionReturn(0);
980316643e7SJed Brown }
981316643e7SJed Brown 
9824a2ae208SSatish Balay #undef __FUNCT__
9834a2ae208SSatish Balay #define __FUNCT__ "TSView"
9847e2c5f70SBarry Smith /*@C
985d763cef2SBarry Smith     TSView - Prints the TS data structure.
986d763cef2SBarry Smith 
9874c49b128SBarry Smith     Collective on TS
988d763cef2SBarry Smith 
989d763cef2SBarry Smith     Input Parameters:
990d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
991d763cef2SBarry Smith -   viewer - visualization context
992d763cef2SBarry Smith 
993d763cef2SBarry Smith     Options Database Key:
994d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
995d763cef2SBarry Smith 
996d763cef2SBarry Smith     Notes:
997d763cef2SBarry Smith     The available visualization contexts include
998b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
999b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1000d763cef2SBarry Smith          output where only the first processor opens
1001d763cef2SBarry Smith          the file.  All other processors send their
1002d763cef2SBarry Smith          data to the first processor to print.
1003d763cef2SBarry Smith 
1004d763cef2SBarry Smith     The user can open an alternative visualization context with
1005b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
1006d763cef2SBarry Smith 
1007d763cef2SBarry Smith     Level: beginner
1008d763cef2SBarry Smith 
1009d763cef2SBarry Smith .keywords: TS, timestep, view
1010d763cef2SBarry Smith 
1011b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
1012d763cef2SBarry Smith @*/
10137087cfbeSBarry Smith PetscErrorCode  TSView(TS ts,PetscViewer viewer)
1014d763cef2SBarry Smith {
1015dfbe8321SBarry Smith   PetscErrorCode ierr;
101619fd82e9SBarry Smith   TSType         type;
1017089b2837SJed Brown   PetscBool      iascii,isstring,isundials;
1018d763cef2SBarry Smith 
1019d763cef2SBarry Smith   PetscFunctionBegin;
10200700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10213050cee2SBarry Smith   if (!viewer) {
10227adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
10233050cee2SBarry Smith   }
10240700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1025c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
1026fd16b177SBarry Smith 
1027251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1028251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
102932077d6dSBarry Smith   if (iascii) {
1030317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr);
103177431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
1032a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
1033d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
10345ef26d82SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);CHKERRQ(ierr);
1035c610991cSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr);
1036d763cef2SBarry Smith     }
10375ef26d82SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);CHKERRQ(ierr);
1038193ac0bcSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
1039d52bd9f3SBarry Smith     if (ts->ops->view) {
1040d52bd9f3SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1041d52bd9f3SBarry Smith       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
1042d52bd9f3SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1043d52bd9f3SBarry Smith     }
10440f5bd95cSBarry Smith   } else if (isstring) {
1045a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
1046b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
1047d763cef2SBarry Smith   }
1048b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1049251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
1050b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1051d763cef2SBarry Smith   PetscFunctionReturn(0);
1052d763cef2SBarry Smith }
1053d763cef2SBarry Smith 
1054d763cef2SBarry Smith 
10554a2ae208SSatish Balay #undef __FUNCT__
10564a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
1057b07ff414SBarry Smith /*@
1058d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
1059d763cef2SBarry Smith    the timesteppers.
1060d763cef2SBarry Smith 
10613f9fe445SBarry Smith    Logically Collective on TS
1062d763cef2SBarry Smith 
1063d763cef2SBarry Smith    Input Parameters:
1064d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1065d763cef2SBarry Smith -  usrP - optional user context
1066d763cef2SBarry Smith 
1067d763cef2SBarry Smith    Level: intermediate
1068d763cef2SBarry Smith 
1069d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
1070d763cef2SBarry Smith 
1071d763cef2SBarry Smith .seealso: TSGetApplicationContext()
1072d763cef2SBarry Smith @*/
10737087cfbeSBarry Smith PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
1074d763cef2SBarry Smith {
1075d763cef2SBarry Smith   PetscFunctionBegin;
10760700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1077d763cef2SBarry Smith   ts->user = usrP;
1078d763cef2SBarry Smith   PetscFunctionReturn(0);
1079d763cef2SBarry Smith }
1080d763cef2SBarry Smith 
10814a2ae208SSatish Balay #undef __FUNCT__
10824a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
1083b07ff414SBarry Smith /*@
1084d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
1085d763cef2SBarry Smith     timestepper.
1086d763cef2SBarry Smith 
1087d763cef2SBarry Smith     Not Collective
1088d763cef2SBarry Smith 
1089d763cef2SBarry Smith     Input Parameter:
1090d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
1091d763cef2SBarry Smith 
1092d763cef2SBarry Smith     Output Parameter:
1093d763cef2SBarry Smith .   usrP - user context
1094d763cef2SBarry Smith 
1095d763cef2SBarry Smith     Level: intermediate
1096d763cef2SBarry Smith 
1097d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
1098d763cef2SBarry Smith 
1099d763cef2SBarry Smith .seealso: TSSetApplicationContext()
1100d763cef2SBarry Smith @*/
1101e71120c6SJed Brown PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
1102d763cef2SBarry Smith {
1103d763cef2SBarry Smith   PetscFunctionBegin;
11040700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1105e71120c6SJed Brown   *(void**)usrP = ts->user;
1106d763cef2SBarry Smith   PetscFunctionReturn(0);
1107d763cef2SBarry Smith }
1108d763cef2SBarry Smith 
11094a2ae208SSatish Balay #undef __FUNCT__
11104a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
1111d763cef2SBarry Smith /*@
1112b8123daeSJed Brown    TSGetTimeStepNumber - Gets the number of time steps completed.
1113d763cef2SBarry Smith 
1114d763cef2SBarry Smith    Not Collective
1115d763cef2SBarry Smith 
1116d763cef2SBarry Smith    Input Parameter:
1117d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1118d763cef2SBarry Smith 
1119d763cef2SBarry Smith    Output Parameter:
1120b8123daeSJed Brown .  iter - number of steps completed so far
1121d763cef2SBarry Smith 
1122d763cef2SBarry Smith    Level: intermediate
1123d763cef2SBarry Smith 
1124d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
1125b8123daeSJed Brown .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStep()
1126d763cef2SBarry Smith @*/
11277087cfbeSBarry Smith PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
1128d763cef2SBarry Smith {
1129d763cef2SBarry Smith   PetscFunctionBegin;
11300700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
11314482741eSBarry Smith   PetscValidIntPointer(iter,2);
1132d763cef2SBarry Smith   *iter = ts->steps;
1133d763cef2SBarry Smith   PetscFunctionReturn(0);
1134d763cef2SBarry Smith }
1135d763cef2SBarry Smith 
11364a2ae208SSatish Balay #undef __FUNCT__
11374a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
1138d763cef2SBarry Smith /*@
1139d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
1140d763cef2SBarry Smith    as well as the initial time.
1141d763cef2SBarry Smith 
11423f9fe445SBarry Smith    Logically Collective on TS
1143d763cef2SBarry Smith 
1144d763cef2SBarry Smith    Input Parameters:
1145d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1146d763cef2SBarry Smith .  initial_time - the initial time
1147d763cef2SBarry Smith -  time_step - the size of the timestep
1148d763cef2SBarry Smith 
1149d763cef2SBarry Smith    Level: intermediate
1150d763cef2SBarry Smith 
1151d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
1152d763cef2SBarry Smith 
1153d763cef2SBarry Smith .keywords: TS, set, initial, timestep
1154d763cef2SBarry Smith @*/
11557087cfbeSBarry Smith PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
1156d763cef2SBarry Smith {
1157e144a568SJed Brown   PetscErrorCode ierr;
1158e144a568SJed Brown 
1159d763cef2SBarry Smith   PetscFunctionBegin;
11600700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1161e144a568SJed Brown   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
1162d8cd7023SBarry Smith   ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr);
1163d763cef2SBarry Smith   PetscFunctionReturn(0);
1164d763cef2SBarry Smith }
1165d763cef2SBarry Smith 
11664a2ae208SSatish Balay #undef __FUNCT__
11674a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
1168d763cef2SBarry Smith /*@
1169d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
1170d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
1171d763cef2SBarry Smith 
11723f9fe445SBarry Smith    Logically Collective on TS
1173d763cef2SBarry Smith 
1174d763cef2SBarry Smith    Input Parameters:
1175d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1176d763cef2SBarry Smith -  time_step - the size of the timestep
1177d763cef2SBarry Smith 
1178d763cef2SBarry Smith    Level: intermediate
1179d763cef2SBarry Smith 
1180d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1181d763cef2SBarry Smith 
1182d763cef2SBarry Smith .keywords: TS, set, timestep
1183d763cef2SBarry Smith @*/
11847087cfbeSBarry Smith PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
1185d763cef2SBarry Smith {
1186d763cef2SBarry Smith   PetscFunctionBegin;
11870700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1188c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,time_step,2);
1189d763cef2SBarry Smith   ts->time_step = time_step;
119031748224SBarry Smith   ts->time_step_orig = time_step;
1191d763cef2SBarry Smith   PetscFunctionReturn(0);
1192d763cef2SBarry Smith }
1193d763cef2SBarry Smith 
11944a2ae208SSatish Balay #undef __FUNCT__
1195a43b19c4SJed Brown #define __FUNCT__ "TSSetExactFinalTime"
1196a43b19c4SJed Brown /*@
1197a43b19c4SJed Brown    TSSetExactFinalTime - Determines whether to interpolate solution to the
1198a43b19c4SJed Brown       exact final time requested by the user or just returns it at the final time
1199a43b19c4SJed Brown       it computed.
1200a43b19c4SJed Brown 
1201a43b19c4SJed Brown   Logically Collective on TS
1202a43b19c4SJed Brown 
1203a43b19c4SJed Brown    Input Parameter:
1204a43b19c4SJed Brown +   ts - the time-step context
1205a43b19c4SJed Brown -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
1206a43b19c4SJed Brown 
1207a43b19c4SJed Brown    Level: beginner
1208a43b19c4SJed Brown 
1209a43b19c4SJed Brown .seealso: TSSetDuration()
1210a43b19c4SJed Brown @*/
1211a43b19c4SJed Brown PetscErrorCode  TSSetExactFinalTime(TS ts,PetscBool flg)
1212a43b19c4SJed Brown {
1213a43b19c4SJed Brown 
1214a43b19c4SJed Brown   PetscFunctionBegin;
1215a43b19c4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1216d8cd7023SBarry Smith   PetscValidLogicalCollectiveBool(ts,flg,2);
1217a43b19c4SJed Brown   ts->exact_final_time = flg;
1218a43b19c4SJed Brown   PetscFunctionReturn(0);
1219a43b19c4SJed Brown }
1220a43b19c4SJed Brown 
1221a43b19c4SJed Brown #undef __FUNCT__
12224a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
1223d763cef2SBarry Smith /*@
1224d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
1225d763cef2SBarry Smith 
1226d763cef2SBarry Smith    Not Collective
1227d763cef2SBarry Smith 
1228d763cef2SBarry Smith    Input Parameter:
1229d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1230d763cef2SBarry Smith 
1231d763cef2SBarry Smith    Output Parameter:
1232d763cef2SBarry Smith .  dt - the current timestep size
1233d763cef2SBarry Smith 
1234d763cef2SBarry Smith    Level: intermediate
1235d763cef2SBarry Smith 
1236d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1237d763cef2SBarry Smith 
1238d763cef2SBarry Smith .keywords: TS, get, timestep
1239d763cef2SBarry Smith @*/
12407087cfbeSBarry Smith PetscErrorCode  TSGetTimeStep(TS ts,PetscReal* dt)
1241d763cef2SBarry Smith {
1242d763cef2SBarry Smith   PetscFunctionBegin;
12430700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1244f7cf8827SBarry Smith   PetscValidRealPointer(dt,2);
1245d763cef2SBarry Smith   *dt = ts->time_step;
1246d763cef2SBarry Smith   PetscFunctionReturn(0);
1247d763cef2SBarry Smith }
1248d763cef2SBarry Smith 
12494a2ae208SSatish Balay #undef __FUNCT__
12504a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
1251d8e5e3e6SSatish Balay /*@
1252d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
1253d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
1254d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
1255d763cef2SBarry Smith    the solution at the next timestep has been calculated.
1256d763cef2SBarry Smith 
1257d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
1258d763cef2SBarry Smith 
1259d763cef2SBarry Smith    Input Parameter:
1260d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1261d763cef2SBarry Smith 
1262d763cef2SBarry Smith    Output Parameter:
1263d763cef2SBarry Smith .  v - the vector containing the solution
1264d763cef2SBarry Smith 
1265d763cef2SBarry Smith    Level: intermediate
1266d763cef2SBarry Smith 
1267d763cef2SBarry Smith .seealso: TSGetTimeStep()
1268d763cef2SBarry Smith 
1269d763cef2SBarry Smith .keywords: TS, timestep, get, solution
1270d763cef2SBarry Smith @*/
12717087cfbeSBarry Smith PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1272d763cef2SBarry Smith {
1273d763cef2SBarry Smith   PetscFunctionBegin;
12740700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12754482741eSBarry Smith   PetscValidPointer(v,2);
12768737fe31SLisandro Dalcin   *v = ts->vec_sol;
1277d763cef2SBarry Smith   PetscFunctionReturn(0);
1278d763cef2SBarry Smith }
1279d763cef2SBarry Smith 
1280bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
12814a2ae208SSatish Balay #undef __FUNCT__
1282bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
1283d8e5e3e6SSatish Balay /*@
1284bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
1285d763cef2SBarry Smith 
1286bdad233fSMatthew Knepley   Not collective
1287d763cef2SBarry Smith 
1288bdad233fSMatthew Knepley   Input Parameters:
1289bdad233fSMatthew Knepley + ts   - The TS
1290bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1291d763cef2SBarry Smith .vb
1292d763cef2SBarry Smith          U_t = A U
1293d763cef2SBarry Smith          U_t = A(t) U
1294d763cef2SBarry Smith          U_t = F(t,U)
1295d763cef2SBarry Smith .ve
1296d763cef2SBarry Smith 
1297d763cef2SBarry Smith    Level: beginner
1298d763cef2SBarry Smith 
1299bdad233fSMatthew Knepley .keywords: TS, problem type
1300bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1301d763cef2SBarry Smith @*/
13027087cfbeSBarry Smith PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1303a7cc72afSBarry Smith {
13049e2a6581SJed Brown   PetscErrorCode ierr;
13059e2a6581SJed Brown 
1306d763cef2SBarry Smith   PetscFunctionBegin;
13070700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1308bdad233fSMatthew Knepley   ts->problem_type = type;
13099e2a6581SJed Brown   if (type == TS_LINEAR) {
13109e2a6581SJed Brown     SNES snes;
13119e2a6581SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
13129e2a6581SJed Brown     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
13139e2a6581SJed Brown   }
1314d763cef2SBarry Smith   PetscFunctionReturn(0);
1315d763cef2SBarry Smith }
1316d763cef2SBarry Smith 
1317bdad233fSMatthew Knepley #undef __FUNCT__
1318bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
1319bdad233fSMatthew Knepley /*@C
1320bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
1321bdad233fSMatthew Knepley 
1322bdad233fSMatthew Knepley   Not collective
1323bdad233fSMatthew Knepley 
1324bdad233fSMatthew Knepley   Input Parameter:
1325bdad233fSMatthew Knepley . ts   - The TS
1326bdad233fSMatthew Knepley 
1327bdad233fSMatthew Knepley   Output Parameter:
1328bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1329bdad233fSMatthew Knepley .vb
1330089b2837SJed Brown          M U_t = A U
1331089b2837SJed Brown          M(t) U_t = A(t) U
1332bdad233fSMatthew Knepley          U_t = F(t,U)
1333bdad233fSMatthew Knepley .ve
1334bdad233fSMatthew Knepley 
1335bdad233fSMatthew Knepley    Level: beginner
1336bdad233fSMatthew Knepley 
1337bdad233fSMatthew Knepley .keywords: TS, problem type
1338bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1339bdad233fSMatthew Knepley @*/
13407087cfbeSBarry Smith PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1341a7cc72afSBarry Smith {
1342bdad233fSMatthew Knepley   PetscFunctionBegin;
13430700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
13444482741eSBarry Smith   PetscValidIntPointer(type,2);
1345bdad233fSMatthew Knepley   *type = ts->problem_type;
1346bdad233fSMatthew Knepley   PetscFunctionReturn(0);
1347bdad233fSMatthew Knepley }
1348d763cef2SBarry Smith 
13494a2ae208SSatish Balay #undef __FUNCT__
13504a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
1351d763cef2SBarry Smith /*@
1352d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
1353d763cef2SBarry Smith    of a timestepper.
1354d763cef2SBarry Smith 
1355d763cef2SBarry Smith    Collective on TS
1356d763cef2SBarry Smith 
1357d763cef2SBarry Smith    Input Parameter:
1358d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1359d763cef2SBarry Smith 
1360d763cef2SBarry Smith    Notes:
1361d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
1362d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
1363d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
1364d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
1365d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
1366d763cef2SBarry Smith 
1367d763cef2SBarry Smith    Level: advanced
1368d763cef2SBarry Smith 
1369d763cef2SBarry Smith .keywords: TS, timestep, setup
1370d763cef2SBarry Smith 
1371d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
1372d763cef2SBarry Smith @*/
13737087cfbeSBarry Smith PetscErrorCode  TSSetUp(TS ts)
1374d763cef2SBarry Smith {
1375dfbe8321SBarry Smith   PetscErrorCode ierr;
13766c6b9e74SPeter Brune   DM             dm;
13776c6b9e74SPeter Brune   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
13786c6b9e74SPeter Brune   PetscErrorCode (*jac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
13796c6b9e74SPeter Brune   TSIJacobian    ijac;
13806c6b9e74SPeter Brune   TSRHSJacobian  rhsjac;
1381d763cef2SBarry Smith 
1382d763cef2SBarry Smith   PetscFunctionBegin;
13830700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1384277b19d0SLisandro Dalcin   if (ts->setupcalled) PetscFunctionReturn(0);
1385277b19d0SLisandro Dalcin 
13867adad957SLisandro Dalcin   if (!((PetscObject)ts)->type_name) {
13879596e0b4SJed Brown     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1388d763cef2SBarry Smith   }
1389a43b19c4SJed Brown   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1390277b19d0SLisandro Dalcin 
1391277b19d0SLisandro Dalcin   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1392277b19d0SLisandro Dalcin 
13931c3436cfSJed Brown   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
13941c3436cfSJed Brown 
1395277b19d0SLisandro Dalcin   if (ts->ops->setup) {
1396000e7ae3SMatthew Knepley     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1397277b19d0SLisandro Dalcin   }
1398277b19d0SLisandro Dalcin 
13996c6b9e74SPeter Brune   /* in the case where we've set a DMTSFunction or what have you, we need the default SNESFunction
14006c6b9e74SPeter Brune    to be set right but can't do it elsewhere due to the overreliance on ctx=ts.
14016c6b9e74SPeter Brune    */
14026c6b9e74SPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
14036c6b9e74SPeter Brune   ierr = DMSNESGetFunction(dm,&func,PETSC_NULL);CHKERRQ(ierr);
14046c6b9e74SPeter Brune   if (!func) {
14056c6b9e74SPeter Brune     ierr =DMSNESSetFunction(dm,SNESTSFormFunction,ts);CHKERRQ(ierr);
14066c6b9e74SPeter Brune   }
14076c6b9e74SPeter 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.
14086c6b9e74SPeter Brune      Otherwise, the SNES will use coloring internally to form the Jacobian.
14096c6b9e74SPeter Brune    */
14106c6b9e74SPeter Brune   ierr = DMSNESGetJacobian(dm,&jac,PETSC_NULL);CHKERRQ(ierr);
14116c6b9e74SPeter Brune   ierr = DMTSGetIJacobian(dm,&ijac,PETSC_NULL);CHKERRQ(ierr);
14126c6b9e74SPeter Brune   ierr = DMTSGetRHSJacobian(dm,&rhsjac,PETSC_NULL);CHKERRQ(ierr);
14136c6b9e74SPeter Brune   if (!jac && (ijac || rhsjac)) {
14146c6b9e74SPeter Brune     ierr = DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);CHKERRQ(ierr);
14156c6b9e74SPeter Brune   }
1416277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_TRUE;
1417277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1418277b19d0SLisandro Dalcin }
1419277b19d0SLisandro Dalcin 
1420277b19d0SLisandro Dalcin #undef __FUNCT__
1421277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset"
1422277b19d0SLisandro Dalcin /*@
1423277b19d0SLisandro Dalcin    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1424277b19d0SLisandro Dalcin 
1425277b19d0SLisandro Dalcin    Collective on TS
1426277b19d0SLisandro Dalcin 
1427277b19d0SLisandro Dalcin    Input Parameter:
1428277b19d0SLisandro Dalcin .  ts - the TS context obtained from TSCreate()
1429277b19d0SLisandro Dalcin 
1430277b19d0SLisandro Dalcin    Level: beginner
1431277b19d0SLisandro Dalcin 
1432277b19d0SLisandro Dalcin .keywords: TS, timestep, reset
1433277b19d0SLisandro Dalcin 
1434277b19d0SLisandro Dalcin .seealso: TSCreate(), TSSetup(), TSDestroy()
1435277b19d0SLisandro Dalcin @*/
1436277b19d0SLisandro Dalcin PetscErrorCode  TSReset(TS ts)
1437277b19d0SLisandro Dalcin {
1438277b19d0SLisandro Dalcin   PetscErrorCode ierr;
1439277b19d0SLisandro Dalcin 
1440277b19d0SLisandro Dalcin   PetscFunctionBegin;
1441277b19d0SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1442277b19d0SLisandro Dalcin   if (ts->ops->reset) {
1443277b19d0SLisandro Dalcin     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1444277b19d0SLisandro Dalcin   }
1445277b19d0SLisandro Dalcin   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
14464e684422SJed Brown   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
14474e684422SJed Brown   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1448214bc6a2SJed Brown   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
14496bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1450e3d84a46SJed Brown   ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
1451e3d84a46SJed Brown   ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
145238637c2eSJed Brown   ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);
1453277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_FALSE;
1454d763cef2SBarry Smith   PetscFunctionReturn(0);
1455d763cef2SBarry Smith }
1456d763cef2SBarry Smith 
14574a2ae208SSatish Balay #undef __FUNCT__
14584a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
1459d8e5e3e6SSatish Balay /*@
1460d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
1461d763cef2SBarry Smith    with TSCreate().
1462d763cef2SBarry Smith 
1463d763cef2SBarry Smith    Collective on TS
1464d763cef2SBarry Smith 
1465d763cef2SBarry Smith    Input Parameter:
1466d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1467d763cef2SBarry Smith 
1468d763cef2SBarry Smith    Level: beginner
1469d763cef2SBarry Smith 
1470d763cef2SBarry Smith .keywords: TS, timestepper, destroy
1471d763cef2SBarry Smith 
1472d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
1473d763cef2SBarry Smith @*/
14746bf464f9SBarry Smith PetscErrorCode  TSDestroy(TS *ts)
1475d763cef2SBarry Smith {
14766849ba73SBarry Smith   PetscErrorCode ierr;
1477d763cef2SBarry Smith 
1478d763cef2SBarry Smith   PetscFunctionBegin;
14796bf464f9SBarry Smith   if (!*ts) PetscFunctionReturn(0);
14806bf464f9SBarry Smith   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
14816bf464f9SBarry Smith   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1482d763cef2SBarry Smith 
14836bf464f9SBarry Smith   ierr = TSReset((*ts));CHKERRQ(ierr);
1484277b19d0SLisandro Dalcin 
1485be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
14866bf464f9SBarry Smith   ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr);
14876bf464f9SBarry Smith   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
14886d4c513bSLisandro Dalcin 
148984df9cb4SJed Brown   ierr = TSAdaptDestroy(&(*ts)->adapt);CHKERRQ(ierr);
14906bf464f9SBarry Smith   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
14916bf464f9SBarry Smith   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
14926bf464f9SBarry Smith   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
14936d4c513bSLisandro Dalcin 
1494a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1495d763cef2SBarry Smith   PetscFunctionReturn(0);
1496d763cef2SBarry Smith }
1497d763cef2SBarry Smith 
14984a2ae208SSatish Balay #undef __FUNCT__
14994a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
1500d8e5e3e6SSatish Balay /*@
1501d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1502d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
1503d763cef2SBarry Smith 
1504d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
1505d763cef2SBarry Smith 
1506d763cef2SBarry Smith    Input Parameter:
1507d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1508d763cef2SBarry Smith 
1509d763cef2SBarry Smith    Output Parameter:
1510d763cef2SBarry Smith .  snes - the nonlinear solver context
1511d763cef2SBarry Smith 
1512d763cef2SBarry Smith    Notes:
1513d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
1514d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
151594b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
1516d763cef2SBarry Smith 
1517d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
1518d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
1519d763cef2SBarry Smith 
1520d763cef2SBarry Smith    Level: beginner
1521d763cef2SBarry Smith 
1522d763cef2SBarry Smith .keywords: timestep, get, SNES
1523d763cef2SBarry Smith @*/
15247087cfbeSBarry Smith PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1525d763cef2SBarry Smith {
1526d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1527d372ba47SLisandro Dalcin 
1528d763cef2SBarry Smith   PetscFunctionBegin;
15290700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
15304482741eSBarry Smith   PetscValidPointer(snes,2);
1531d372ba47SLisandro Dalcin   if (!ts->snes) {
1532d372ba47SLisandro Dalcin     ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
15336c6b9e74SPeter Brune     ierr = SNESSetFunction(ts->snes,PETSC_NULL,SNESTSFormFunction,ts);CHKERRQ(ierr);
1534d372ba47SLisandro Dalcin     ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr);
1535d372ba47SLisandro Dalcin     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
1536496e6a7aSJed Brown     if (ts->dm) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
15379e2a6581SJed Brown     if (ts->problem_type == TS_LINEAR) {
15389e2a6581SJed Brown       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
15399e2a6581SJed Brown     }
1540d372ba47SLisandro Dalcin   }
1541d763cef2SBarry Smith   *snes = ts->snes;
1542d763cef2SBarry Smith   PetscFunctionReturn(0);
1543d763cef2SBarry Smith }
1544d763cef2SBarry Smith 
15454a2ae208SSatish Balay #undef __FUNCT__
154694b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
1547d8e5e3e6SSatish Balay /*@
154894b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
1549d763cef2SBarry Smith    a TS (timestepper) context.
1550d763cef2SBarry Smith 
155194b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1552d763cef2SBarry Smith 
1553d763cef2SBarry Smith    Input Parameter:
1554d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1555d763cef2SBarry Smith 
1556d763cef2SBarry Smith    Output Parameter:
155794b7f48cSBarry Smith .  ksp - the nonlinear solver context
1558d763cef2SBarry Smith 
1559d763cef2SBarry Smith    Notes:
156094b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1561d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1562d763cef2SBarry Smith    KSP and PC contexts as well.
1563d763cef2SBarry Smith 
156494b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
156594b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1566d763cef2SBarry Smith 
1567d763cef2SBarry Smith    Level: beginner
1568d763cef2SBarry Smith 
156994b7f48cSBarry Smith .keywords: timestep, get, KSP
1570d763cef2SBarry Smith @*/
15717087cfbeSBarry Smith PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1572d763cef2SBarry Smith {
1573d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1574089b2837SJed Brown   SNES           snes;
1575d372ba47SLisandro Dalcin 
1576d763cef2SBarry Smith   PetscFunctionBegin;
15770700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
15784482741eSBarry Smith   PetscValidPointer(ksp,2);
157917186662SBarry Smith   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1580e32f2f54SBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1581089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1582089b2837SJed Brown   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
1583d763cef2SBarry Smith   PetscFunctionReturn(0);
1584d763cef2SBarry Smith }
1585d763cef2SBarry Smith 
1586d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1587d763cef2SBarry Smith 
15884a2ae208SSatish Balay #undef __FUNCT__
1589adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1590adb62b0dSMatthew Knepley /*@
1591adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1592adb62b0dSMatthew Knepley    maximum time for iteration.
1593adb62b0dSMatthew Knepley 
15943f9fe445SBarry Smith    Not Collective
1595adb62b0dSMatthew Knepley 
1596adb62b0dSMatthew Knepley    Input Parameters:
1597adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1598adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1599adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1600adb62b0dSMatthew Knepley 
1601adb62b0dSMatthew Knepley    Level: intermediate
1602adb62b0dSMatthew Knepley 
1603adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1604adb62b0dSMatthew Knepley @*/
16057087cfbeSBarry Smith PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1606adb62b0dSMatthew Knepley {
1607adb62b0dSMatthew Knepley   PetscFunctionBegin;
16080700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1609abc0a331SBarry Smith   if (maxsteps) {
16104482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1611adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1612adb62b0dSMatthew Knepley   }
1613abc0a331SBarry Smith   if (maxtime) {
16144482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1615adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1616adb62b0dSMatthew Knepley   }
1617adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1618adb62b0dSMatthew Knepley }
1619adb62b0dSMatthew Knepley 
1620adb62b0dSMatthew Knepley #undef __FUNCT__
16214a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1622d763cef2SBarry Smith /*@
1623d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1624d763cef2SBarry Smith    maximum time for iteration.
1625d763cef2SBarry Smith 
16263f9fe445SBarry Smith    Logically Collective on TS
1627d763cef2SBarry Smith 
1628d763cef2SBarry Smith    Input Parameters:
1629d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1630d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1631d763cef2SBarry Smith -  maxtime - final time to iterate to
1632d763cef2SBarry Smith 
1633d763cef2SBarry Smith    Options Database Keys:
1634d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
16353bca7d26SBarry Smith .  -ts_final_time <maxtime> - Sets maxtime
1636d763cef2SBarry Smith 
1637d763cef2SBarry Smith    Notes:
1638d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1639d763cef2SBarry Smith 
1640d763cef2SBarry Smith    Level: intermediate
1641d763cef2SBarry Smith 
1642d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1643a43b19c4SJed Brown 
1644a43b19c4SJed Brown .seealso: TSSetExactFinalTime()
1645d763cef2SBarry Smith @*/
16467087cfbeSBarry Smith PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1647d763cef2SBarry Smith {
1648d763cef2SBarry Smith   PetscFunctionBegin;
16490700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1650c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1651c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,maxtime,2);
165239b7ec4bSSean Farley   if (maxsteps >= 0) ts->max_steps = maxsteps;
165339b7ec4bSSean Farley   if (maxtime != PETSC_DEFAULT) ts->max_time  = maxtime;
1654d763cef2SBarry Smith   PetscFunctionReturn(0);
1655d763cef2SBarry Smith }
1656d763cef2SBarry Smith 
16574a2ae208SSatish Balay #undef __FUNCT__
16584a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1659d763cef2SBarry Smith /*@
1660d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1661d763cef2SBarry Smith    for use by the TS routines.
1662d763cef2SBarry Smith 
16633f9fe445SBarry Smith    Logically Collective on TS and Vec
1664d763cef2SBarry Smith 
1665d763cef2SBarry Smith    Input Parameters:
1666d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1667d763cef2SBarry Smith -  x - the solution vector
1668d763cef2SBarry Smith 
1669d763cef2SBarry Smith    Level: beginner
1670d763cef2SBarry Smith 
1671d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1672d763cef2SBarry Smith @*/
16737087cfbeSBarry Smith PetscErrorCode  TSSetSolution(TS ts,Vec x)
1674d763cef2SBarry Smith {
16758737fe31SLisandro Dalcin   PetscErrorCode ierr;
1676496e6a7aSJed Brown   DM             dm;
16778737fe31SLisandro Dalcin 
1678d763cef2SBarry Smith   PetscFunctionBegin;
16790700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
16800700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
16818737fe31SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
16826bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
16838737fe31SLisandro Dalcin   ts->vec_sol = x;
1684496e6a7aSJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1685496e6a7aSJed Brown   ierr = DMShellSetGlobalVector(dm,x);CHKERRQ(ierr);
1686d763cef2SBarry Smith   PetscFunctionReturn(0);
1687d763cef2SBarry Smith }
1688d763cef2SBarry Smith 
1689e74ef692SMatthew Knepley #undef __FUNCT__
1690e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1691ac226902SBarry Smith /*@C
1692000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
16933f2090d5SJed Brown   called once at the beginning of each time step.
1694000e7ae3SMatthew Knepley 
16953f9fe445SBarry Smith   Logically Collective on TS
1696000e7ae3SMatthew Knepley 
1697000e7ae3SMatthew Knepley   Input Parameters:
1698000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1699000e7ae3SMatthew Knepley - func - The function
1700000e7ae3SMatthew Knepley 
1701000e7ae3SMatthew Knepley   Calling sequence of func:
1702000e7ae3SMatthew Knepley . func (TS ts);
1703000e7ae3SMatthew Knepley 
1704000e7ae3SMatthew Knepley   Level: intermediate
1705000e7ae3SMatthew Knepley 
1706b8123daeSJed Brown   Note:
1707b8123daeSJed Brown   If a step is rejected, TSStep() will call this routine again before each attempt.
1708b8123daeSJed Brown   The last completed time step number can be queried using TSGetTimeStepNumber(), the
1709b8123daeSJed Brown   size of the step being attempted can be obtained using TSGetTimeStep().
1710b8123daeSJed Brown 
1711000e7ae3SMatthew Knepley .keywords: TS, timestep
1712b8123daeSJed Brown .seealso: TSSetPreStage(), TSSetPostStep(), TSStep()
1713000e7ae3SMatthew Knepley @*/
17147087cfbeSBarry Smith PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1715000e7ae3SMatthew Knepley {
1716000e7ae3SMatthew Knepley   PetscFunctionBegin;
17170700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1718000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1719000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1720000e7ae3SMatthew Knepley }
1721000e7ae3SMatthew Knepley 
1722e74ef692SMatthew Knepley #undef __FUNCT__
17233f2090d5SJed Brown #define __FUNCT__ "TSPreStep"
172409ee8438SJed Brown /*@
17253f2090d5SJed Brown   TSPreStep - Runs the user-defined pre-step function.
17263f2090d5SJed Brown 
17273f2090d5SJed Brown   Collective on TS
17283f2090d5SJed Brown 
17293f2090d5SJed Brown   Input Parameters:
17303f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
17313f2090d5SJed Brown 
17323f2090d5SJed Brown   Notes:
17333f2090d5SJed Brown   TSPreStep() is typically used within time stepping implementations,
17343f2090d5SJed Brown   so most users would not generally call this routine themselves.
17353f2090d5SJed Brown 
17363f2090d5SJed Brown   Level: developer
17373f2090d5SJed Brown 
17383f2090d5SJed Brown .keywords: TS, timestep
1739b8123daeSJed Brown .seealso: TSSetPreStep(), TSPreStage(), TSPostStep()
17403f2090d5SJed Brown @*/
17417087cfbeSBarry Smith PetscErrorCode  TSPreStep(TS ts)
17423f2090d5SJed Brown {
17433f2090d5SJed Brown   PetscErrorCode ierr;
17443f2090d5SJed Brown 
17453f2090d5SJed Brown   PetscFunctionBegin;
17460700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
174772ac3e02SJed Brown   if (ts->ops->prestep) {
17483f2090d5SJed Brown     PetscStackPush("TS PreStep function");
17493f2090d5SJed Brown     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
17503f2090d5SJed Brown     PetscStackPop;
1751312ce896SJed Brown   }
17523f2090d5SJed Brown   PetscFunctionReturn(0);
17533f2090d5SJed Brown }
17543f2090d5SJed Brown 
17553f2090d5SJed Brown #undef __FUNCT__
1756b8123daeSJed Brown #define __FUNCT__ "TSSetPreStage"
1757b8123daeSJed Brown /*@C
1758b8123daeSJed Brown   TSSetPreStage - Sets the general-purpose function
1759b8123daeSJed Brown   called once at the beginning of each stage.
1760b8123daeSJed Brown 
1761b8123daeSJed Brown   Logically Collective on TS
1762b8123daeSJed Brown 
1763b8123daeSJed Brown   Input Parameters:
1764b8123daeSJed Brown + ts   - The TS context obtained from TSCreate()
1765b8123daeSJed Brown - func - The function
1766b8123daeSJed Brown 
1767b8123daeSJed Brown   Calling sequence of func:
1768b8123daeSJed Brown . PetscErrorCode func(TS ts, PetscReal stagetime);
1769b8123daeSJed Brown 
1770b8123daeSJed Brown   Level: intermediate
1771b8123daeSJed Brown 
1772b8123daeSJed Brown   Note:
1773b8123daeSJed Brown   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
1774b8123daeSJed Brown   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
1775b8123daeSJed Brown   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
1776b8123daeSJed Brown 
1777b8123daeSJed Brown .keywords: TS, timestep
1778b8123daeSJed Brown .seealso: TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
1779b8123daeSJed Brown @*/
1780b8123daeSJed Brown PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
1781b8123daeSJed Brown {
1782b8123daeSJed Brown   PetscFunctionBegin;
1783b8123daeSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1784b8123daeSJed Brown   ts->ops->prestage = func;
1785b8123daeSJed Brown   PetscFunctionReturn(0);
1786b8123daeSJed Brown }
1787b8123daeSJed Brown 
1788b8123daeSJed Brown #undef __FUNCT__
1789b8123daeSJed Brown #define __FUNCT__ "TSPreStage"
1790b8123daeSJed Brown /*@
1791b8123daeSJed Brown   TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()
1792b8123daeSJed Brown 
1793b8123daeSJed Brown   Collective on TS
1794b8123daeSJed Brown 
1795b8123daeSJed Brown   Input Parameters:
1796b8123daeSJed Brown . ts   - The TS context obtained from TSCreate()
1797b8123daeSJed Brown 
1798b8123daeSJed Brown   Notes:
1799b8123daeSJed Brown   TSPreStage() is typically used within time stepping implementations,
1800b8123daeSJed Brown   most users would not generally call this routine themselves.
1801b8123daeSJed Brown 
1802b8123daeSJed Brown   Level: developer
1803b8123daeSJed Brown 
1804b8123daeSJed Brown .keywords: TS, timestep
1805b8123daeSJed Brown .seealso: TSSetPreStep(), TSPreStep(), TSPostStep()
1806b8123daeSJed Brown @*/
1807b8123daeSJed Brown PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
1808b8123daeSJed Brown {
1809b8123daeSJed Brown   PetscErrorCode ierr;
1810b8123daeSJed Brown 
1811b8123daeSJed Brown   PetscFunctionBegin;
1812b8123daeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1813b8123daeSJed Brown   if (ts->ops->prestage) {
1814b8123daeSJed Brown     PetscStackPush("TS PreStage function");
1815b8123daeSJed Brown     ierr = (*ts->ops->prestage)(ts,stagetime);CHKERRQ(ierr);
1816b8123daeSJed Brown     PetscStackPop;
1817b8123daeSJed Brown   }
1818b8123daeSJed Brown   PetscFunctionReturn(0);
1819b8123daeSJed Brown }
1820b8123daeSJed Brown 
1821b8123daeSJed Brown #undef __FUNCT__
1822e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1823ac226902SBarry Smith /*@C
1824000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
18253f2090d5SJed Brown   called once at the end of each time step.
1826000e7ae3SMatthew Knepley 
18273f9fe445SBarry Smith   Logically Collective on TS
1828000e7ae3SMatthew Knepley 
1829000e7ae3SMatthew Knepley   Input Parameters:
1830000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1831000e7ae3SMatthew Knepley - func - The function
1832000e7ae3SMatthew Knepley 
1833000e7ae3SMatthew Knepley   Calling sequence of func:
1834b8123daeSJed Brown $ func (TS ts);
1835000e7ae3SMatthew Knepley 
1836000e7ae3SMatthew Knepley   Level: intermediate
1837000e7ae3SMatthew Knepley 
1838000e7ae3SMatthew Knepley .keywords: TS, timestep
1839b8123daeSJed Brown .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime()
1840000e7ae3SMatthew Knepley @*/
18417087cfbeSBarry Smith PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1842000e7ae3SMatthew Knepley {
1843000e7ae3SMatthew Knepley   PetscFunctionBegin;
18440700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1845000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1846000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1847000e7ae3SMatthew Knepley }
1848000e7ae3SMatthew Knepley 
1849e74ef692SMatthew Knepley #undef __FUNCT__
18503f2090d5SJed Brown #define __FUNCT__ "TSPostStep"
185109ee8438SJed Brown /*@
18523f2090d5SJed Brown   TSPostStep - Runs the user-defined post-step function.
18533f2090d5SJed Brown 
18543f2090d5SJed Brown   Collective on TS
18553f2090d5SJed Brown 
18563f2090d5SJed Brown   Input Parameters:
18573f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
18583f2090d5SJed Brown 
18593f2090d5SJed Brown   Notes:
18603f2090d5SJed Brown   TSPostStep() is typically used within time stepping implementations,
18613f2090d5SJed Brown   so most users would not generally call this routine themselves.
18623f2090d5SJed Brown 
18633f2090d5SJed Brown   Level: developer
18643f2090d5SJed Brown 
18653f2090d5SJed Brown .keywords: TS, timestep
18663f2090d5SJed Brown @*/
18677087cfbeSBarry Smith PetscErrorCode  TSPostStep(TS ts)
18683f2090d5SJed Brown {
18693f2090d5SJed Brown   PetscErrorCode ierr;
18703f2090d5SJed Brown 
18713f2090d5SJed Brown   PetscFunctionBegin;
18720700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
187372ac3e02SJed Brown   if (ts->ops->poststep) {
18743f2090d5SJed Brown     PetscStackPush("TS PostStep function");
18753f2090d5SJed Brown     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
18763f2090d5SJed Brown     PetscStackPop;
187772ac3e02SJed Brown   }
18783f2090d5SJed Brown   PetscFunctionReturn(0);
18793f2090d5SJed Brown }
18803f2090d5SJed Brown 
1881d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1882d763cef2SBarry Smith 
18834a2ae208SSatish Balay #undef __FUNCT__
1884a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet"
1885d763cef2SBarry Smith /*@C
1886a6570f20SBarry Smith    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1887d763cef2SBarry Smith    timestep to display the iteration's  progress.
1888d763cef2SBarry Smith 
18893f9fe445SBarry Smith    Logically Collective on TS
1890d763cef2SBarry Smith 
1891d763cef2SBarry Smith    Input Parameters:
1892d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1893e213d8f1SJed Brown .  monitor - monitoring routine
1894329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1895b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1896b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1897b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1898d763cef2SBarry Smith 
1899e213d8f1SJed Brown    Calling sequence of monitor:
1900e213d8f1SJed Brown $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1901d763cef2SBarry Smith 
1902d763cef2SBarry Smith +    ts - the TS context
1903d763cef2SBarry Smith .    steps - iteration number
19041f06c33eSBarry Smith .    time - current time
1905d763cef2SBarry Smith .    x - current iterate
1906d763cef2SBarry Smith -    mctx - [optional] monitoring context
1907d763cef2SBarry Smith 
1908d763cef2SBarry Smith    Notes:
1909d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1910d763cef2SBarry Smith    already has been loaded.
1911d763cef2SBarry Smith 
1912025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each TS object
1913025f1a04SBarry Smith 
1914d763cef2SBarry Smith    Level: intermediate
1915d763cef2SBarry Smith 
1916d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1917d763cef2SBarry Smith 
1918a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel()
1919d763cef2SBarry Smith @*/
1920c2efdce3SBarry Smith PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1921d763cef2SBarry Smith {
1922d763cef2SBarry Smith   PetscFunctionBegin;
19230700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
192417186662SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1925d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1926329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1927d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1928d763cef2SBarry Smith   PetscFunctionReturn(0);
1929d763cef2SBarry Smith }
1930d763cef2SBarry Smith 
19314a2ae208SSatish Balay #undef __FUNCT__
1932a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel"
1933d763cef2SBarry Smith /*@C
1934a6570f20SBarry Smith    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1935d763cef2SBarry Smith 
19363f9fe445SBarry Smith    Logically Collective on TS
1937d763cef2SBarry Smith 
1938d763cef2SBarry Smith    Input Parameters:
1939d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1940d763cef2SBarry Smith 
1941d763cef2SBarry Smith    Notes:
1942d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1943d763cef2SBarry Smith 
1944d763cef2SBarry Smith    Level: intermediate
1945d763cef2SBarry Smith 
1946d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1947d763cef2SBarry Smith 
1948a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet()
1949d763cef2SBarry Smith @*/
19507087cfbeSBarry Smith PetscErrorCode  TSMonitorCancel(TS ts)
1951d763cef2SBarry Smith {
1952d952e501SBarry Smith   PetscErrorCode ierr;
1953d952e501SBarry Smith   PetscInt       i;
1954d952e501SBarry Smith 
1955d763cef2SBarry Smith   PetscFunctionBegin;
19560700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1957d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1958d952e501SBarry Smith     if (ts->mdestroy[i]) {
19593c4aec1bSBarry Smith       ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
1960d952e501SBarry Smith     }
1961d952e501SBarry Smith   }
1962d763cef2SBarry Smith   ts->numbermonitors = 0;
1963d763cef2SBarry Smith   PetscFunctionReturn(0);
1964d763cef2SBarry Smith }
1965d763cef2SBarry Smith 
19664a2ae208SSatish Balay #undef __FUNCT__
1967a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault"
1968d8e5e3e6SSatish Balay /*@
1969a6570f20SBarry Smith    TSMonitorDefault - Sets the Default monitor
19705516499fSSatish Balay 
19715516499fSSatish Balay    Level: intermediate
197241251cbbSSatish Balay 
19735516499fSSatish Balay .keywords: TS, set, monitor
19745516499fSSatish Balay 
197541251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet()
197641251cbbSSatish Balay @*/
1977649052a6SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1978d763cef2SBarry Smith {
1979dfbe8321SBarry Smith   PetscErrorCode ierr;
1980649052a6SBarry Smith   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1981d132466eSBarry Smith 
1982d763cef2SBarry Smith   PetscFunctionBegin;
1983649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1984971832b6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr);
1985649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1986d763cef2SBarry Smith   PetscFunctionReturn(0);
1987d763cef2SBarry Smith }
1988d763cef2SBarry Smith 
19894a2ae208SSatish Balay #undef __FUNCT__
1990cd652676SJed Brown #define __FUNCT__ "TSSetRetainStages"
1991cd652676SJed Brown /*@
1992cd652676SJed Brown    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
1993cd652676SJed Brown 
1994cd652676SJed Brown    Logically Collective on TS
1995cd652676SJed Brown 
1996cd652676SJed Brown    Input Argument:
1997cd652676SJed Brown .  ts - time stepping context
1998cd652676SJed Brown 
1999cd652676SJed Brown    Output Argument:
2000cd652676SJed Brown .  flg - PETSC_TRUE or PETSC_FALSE
2001cd652676SJed Brown 
2002cd652676SJed Brown    Level: intermediate
2003cd652676SJed Brown 
2004cd652676SJed Brown .keywords: TS, set
2005cd652676SJed Brown 
2006cd652676SJed Brown .seealso: TSInterpolate(), TSSetPostStep()
2007cd652676SJed Brown @*/
2008cd652676SJed Brown PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
2009cd652676SJed Brown {
2010cd652676SJed Brown 
2011cd652676SJed Brown   PetscFunctionBegin;
2012cd652676SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2013cd652676SJed Brown   ts->retain_stages = flg;
2014cd652676SJed Brown   PetscFunctionReturn(0);
2015cd652676SJed Brown }
2016cd652676SJed Brown 
2017cd652676SJed Brown #undef __FUNCT__
2018cd652676SJed Brown #define __FUNCT__ "TSInterpolate"
2019cd652676SJed Brown /*@
2020cd652676SJed Brown    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
2021cd652676SJed Brown 
2022cd652676SJed Brown    Collective on TS
2023cd652676SJed Brown 
2024cd652676SJed Brown    Input Argument:
2025cd652676SJed Brown +  ts - time stepping context
2026cd652676SJed Brown -  t - time to interpolate to
2027cd652676SJed Brown 
2028cd652676SJed Brown    Output Argument:
2029cd652676SJed Brown .  X - state at given time
2030cd652676SJed Brown 
2031cd652676SJed Brown    Notes:
2032cd652676SJed Brown    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
2033cd652676SJed Brown 
2034cd652676SJed Brown    Level: intermediate
2035cd652676SJed Brown 
2036cd652676SJed Brown    Developer Notes:
2037cd652676SJed Brown    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
2038cd652676SJed Brown 
2039cd652676SJed Brown .keywords: TS, set
2040cd652676SJed Brown 
2041cd652676SJed Brown .seealso: TSSetRetainStages(), TSSetPostStep()
2042cd652676SJed Brown @*/
2043cd652676SJed Brown PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
2044cd652676SJed Brown {
2045cd652676SJed Brown   PetscErrorCode ierr;
2046cd652676SJed Brown 
2047cd652676SJed Brown   PetscFunctionBegin;
2048cd652676SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2049d8cd7023SBarry 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);
2050cd652676SJed Brown   if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
2051cd652676SJed Brown   ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr);
2052cd652676SJed Brown   PetscFunctionReturn(0);
2053cd652676SJed Brown }
2054cd652676SJed Brown 
2055cd652676SJed Brown #undef __FUNCT__
20564a2ae208SSatish Balay #define __FUNCT__ "TSStep"
2057d763cef2SBarry Smith /*@
20586d9e5789SSean Farley    TSStep - Steps one time step
2059d763cef2SBarry Smith 
2060d763cef2SBarry Smith    Collective on TS
2061d763cef2SBarry Smith 
2062d763cef2SBarry Smith    Input Parameter:
2063d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
2064d763cef2SBarry Smith 
20656d9e5789SSean Farley    Level: intermediate
2066d763cef2SBarry Smith 
2067b8123daeSJed Brown    Notes:
2068b8123daeSJed Brown    The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
2069b8123daeSJed Brown    be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.
2070b8123daeSJed Brown 
2071d763cef2SBarry Smith .keywords: TS, timestep, solve
2072d763cef2SBarry Smith 
2073b8123daeSJed Brown .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage()
2074d763cef2SBarry Smith @*/
2075193ac0bcSJed Brown PetscErrorCode  TSStep(TS ts)
2076d763cef2SBarry Smith {
2077362cd11cSLisandro Dalcin   PetscReal      ptime_prev;
2078dfbe8321SBarry Smith   PetscErrorCode ierr;
2079d763cef2SBarry Smith 
2080d763cef2SBarry Smith   PetscFunctionBegin;
20810700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
2082d405a339SMatthew Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
2083d405a339SMatthew Knepley 
2084362cd11cSLisandro Dalcin   ts->reason = TS_CONVERGED_ITERATING;
2085362cd11cSLisandro Dalcin 
2086362cd11cSLisandro Dalcin   ptime_prev = ts->ptime;
2087d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
2088193ac0bcSJed Brown   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
2089d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
2090362cd11cSLisandro Dalcin   ts->time_step_prev = ts->ptime - ptime_prev;
2091362cd11cSLisandro Dalcin 
2092362cd11cSLisandro Dalcin   if (ts->reason < 0) {
2093cef5090cSJed Brown     if (ts->errorifstepfailed) {
2094cef5090cSJed Brown       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
2095cef5090cSJed 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]);
2096cef5090cSJed Brown       } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
2097cef5090cSJed Brown     }
2098362cd11cSLisandro Dalcin   } else if (!ts->reason) {
2099362cd11cSLisandro Dalcin     if (ts->steps >= ts->max_steps)
2100362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_ITS;
2101362cd11cSLisandro Dalcin     else if (ts->ptime >= ts->max_time)
2102362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_TIME;
2103362cd11cSLisandro Dalcin   }
2104362cd11cSLisandro Dalcin 
2105d763cef2SBarry Smith   PetscFunctionReturn(0);
2106d763cef2SBarry Smith }
2107d763cef2SBarry Smith 
21084a2ae208SSatish Balay #undef __FUNCT__
210905175c85SJed Brown #define __FUNCT__ "TSEvaluateStep"
211005175c85SJed Brown /*@
211105175c85SJed Brown    TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
211205175c85SJed Brown 
21131c3436cfSJed Brown    Collective on TS
211405175c85SJed Brown 
211505175c85SJed Brown    Input Arguments:
21161c3436cfSJed Brown +  ts - time stepping context
21171c3436cfSJed Brown .  order - desired order of accuracy
21181c3436cfSJed Brown -  done - whether the step was evaluated at this order (pass PETSC_NULL to generate an error if not available)
211905175c85SJed Brown 
212005175c85SJed Brown    Output Arguments:
21211c3436cfSJed Brown .  X - state at the end of the current step
212205175c85SJed Brown 
212305175c85SJed Brown    Level: advanced
212405175c85SJed Brown 
2125108c343cSJed Brown    Notes:
2126108c343cSJed Brown    This function cannot be called until all stages have been evaluated.
2127108c343cSJed 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.
2128108c343cSJed Brown 
21291c3436cfSJed Brown .seealso: TSStep(), TSAdapt
213005175c85SJed Brown @*/
21311c3436cfSJed Brown PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done)
213205175c85SJed Brown {
213305175c85SJed Brown   PetscErrorCode ierr;
213405175c85SJed Brown 
213505175c85SJed Brown   PetscFunctionBegin;
213605175c85SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
213705175c85SJed Brown   PetscValidType(ts,1);
213805175c85SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
21391c3436cfSJed Brown   if (!ts->ops->evaluatestep) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
21401c3436cfSJed Brown   ierr = (*ts->ops->evaluatestep)(ts,order,X,done);CHKERRQ(ierr);
214105175c85SJed Brown   PetscFunctionReturn(0);
214205175c85SJed Brown }
214305175c85SJed Brown 
214405175c85SJed Brown #undef __FUNCT__
21456a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve"
21466a4d4014SLisandro Dalcin /*@
21476a4d4014SLisandro Dalcin    TSSolve - Steps the requested number of timesteps.
21486a4d4014SLisandro Dalcin 
21496a4d4014SLisandro Dalcin    Collective on TS
21506a4d4014SLisandro Dalcin 
21516a4d4014SLisandro Dalcin    Input Parameter:
21526a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
2153362cd11cSLisandro Dalcin -  x - the solution vector
21546a4d4014SLisandro Dalcin 
21555a3a76d0SJed Brown    Output Parameter:
21565a3a76d0SJed Brown .  ftime - time of the state vector x upon completion
21575a3a76d0SJed Brown 
21586a4d4014SLisandro Dalcin    Level: beginner
21596a4d4014SLisandro Dalcin 
21605a3a76d0SJed Brown    Notes:
21615a3a76d0SJed Brown    The final time returned by this function may be different from the time of the internally
21625a3a76d0SJed Brown    held state accessible by TSGetSolution() and TSGetTime() because the method may have
21635a3a76d0SJed Brown    stepped over the final time.
21645a3a76d0SJed Brown 
21656a4d4014SLisandro Dalcin .keywords: TS, timestep, solve
21666a4d4014SLisandro Dalcin 
21676a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep()
21686a4d4014SLisandro Dalcin @*/
21695a3a76d0SJed Brown PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime)
21706a4d4014SLisandro Dalcin {
21714d7d938eSLisandro Dalcin   PetscBool      flg;
21724d7d938eSLisandro Dalcin   char           filename[PETSC_MAX_PATH_LEN];
21734d7d938eSLisandro Dalcin   PetscViewer    viewer;
21746a4d4014SLisandro Dalcin   PetscErrorCode ierr;
2175f22f69f0SBarry Smith 
21766a4d4014SLisandro Dalcin   PetscFunctionBegin;
21770700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2178193ac0bcSJed Brown   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
21795a3a76d0SJed Brown   if (ts->exact_final_time) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
21805a3a76d0SJed Brown     if (!ts->vec_sol || x == ts->vec_sol) {
21815a3a76d0SJed Brown       Vec y;
21825a3a76d0SJed Brown       ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
21835a3a76d0SJed Brown       ierr = TSSetSolution(ts,y);CHKERRQ(ierr);
21845a3a76d0SJed Brown       ierr = VecDestroy(&y);CHKERRQ(ierr); /* grant ownership */
21855a3a76d0SJed Brown     }
2186362cd11cSLisandro Dalcin     ierr = VecCopy(x,ts->vec_sol);CHKERRQ(ierr);
21875a3a76d0SJed Brown   } else {
2188193ac0bcSJed Brown     ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
21895a3a76d0SJed Brown   }
2190b5d403baSSean Farley   ierr = TSSetUp(ts);CHKERRQ(ierr);
21916a4d4014SLisandro Dalcin   /* reset time step and iteration counters */
2192193ac0bcSJed Brown   ts->steps = 0;
21935ef26d82SJed Brown   ts->ksp_its = 0;
21945ef26d82SJed Brown   ts->snes_its = 0;
2195c610991cSLisandro Dalcin   ts->num_snes_failures = 0;
2196c610991cSLisandro Dalcin   ts->reject = 0;
2197193ac0bcSJed Brown   ts->reason = TS_CONVERGED_ITERATING;
2198193ac0bcSJed Brown 
2199193ac0bcSJed Brown   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
2200193ac0bcSJed Brown     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
22015a3a76d0SJed Brown     ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
2202a5b82c69SSean Farley     if (ftime) *ftime = ts->ptime;
2203193ac0bcSJed Brown   } else {
22046a4d4014SLisandro Dalcin     /* steps the requested number of timesteps. */
2205362cd11cSLisandro Dalcin     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2206362cd11cSLisandro Dalcin     if (ts->steps >= ts->max_steps)
2207362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_ITS;
2208362cd11cSLisandro Dalcin     else if (ts->ptime >= ts->max_time)
2209362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_TIME;
2210e1a7a14fSJed Brown     while (!ts->reason) {
2211193ac0bcSJed Brown       ierr = TSStep(ts);CHKERRQ(ierr);
2212193ac0bcSJed Brown       ierr = TSPostStep(ts);CHKERRQ(ierr);
2213193ac0bcSJed Brown       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2214193ac0bcSJed Brown     }
2215362cd11cSLisandro Dalcin     if (ts->exact_final_time && ts->ptime > ts->max_time) {
2216a43b19c4SJed Brown       ierr = TSInterpolate(ts,ts->max_time,x);CHKERRQ(ierr);
22175a3a76d0SJed Brown       if (ftime) *ftime = ts->max_time;
22180574a7fbSJed Brown     } else {
22190574a7fbSJed Brown       ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
22200574a7fbSJed Brown       if (ftime) *ftime = ts->ptime;
22210574a7fbSJed Brown     }
2222193ac0bcSJed Brown   }
22234d7d938eSLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
22244d7d938eSLisandro Dalcin   if (flg && !PetscPreLoadingOn) {
22254d7d938eSLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr);
22264d7d938eSLisandro Dalcin     ierr = TSView(ts,viewer);CHKERRQ(ierr);
22274d7d938eSLisandro Dalcin     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2228193ac0bcSJed Brown   }
22296a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
22306a4d4014SLisandro Dalcin }
22316a4d4014SLisandro Dalcin 
22326a4d4014SLisandro Dalcin #undef __FUNCT__
22334a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
2234228d79bcSJed Brown /*@
2235228d79bcSJed Brown    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
2236228d79bcSJed Brown 
2237228d79bcSJed Brown    Collective on TS
2238228d79bcSJed Brown 
2239228d79bcSJed Brown    Input Parameters:
2240228d79bcSJed Brown +  ts - time stepping context obtained from TSCreate()
2241228d79bcSJed Brown .  step - step number that has just completed
2242228d79bcSJed Brown .  ptime - model time of the state
2243228d79bcSJed Brown -  x - state at the current model time
2244228d79bcSJed Brown 
2245228d79bcSJed Brown    Notes:
2246228d79bcSJed Brown    TSMonitor() is typically used within the time stepping implementations.
2247228d79bcSJed Brown    Users might call this function when using the TSStep() interface instead of TSSolve().
2248228d79bcSJed Brown 
2249228d79bcSJed Brown    Level: advanced
2250228d79bcSJed Brown 
2251228d79bcSJed Brown .keywords: TS, timestep
2252228d79bcSJed Brown @*/
2253a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
2254d763cef2SBarry Smith {
22556849ba73SBarry Smith   PetscErrorCode ierr;
2256a7cc72afSBarry Smith   PetscInt       i,n = ts->numbermonitors;
2257d763cef2SBarry Smith 
2258d763cef2SBarry Smith   PetscFunctionBegin;
2259d763cef2SBarry Smith   for (i=0; i<n; i++) {
2260142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
2261d763cef2SBarry Smith   }
2262d763cef2SBarry Smith   PetscFunctionReturn(0);
2263d763cef2SBarry Smith }
2264d763cef2SBarry Smith 
2265d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
2266d763cef2SBarry Smith 
22674a2ae208SSatish Balay #undef __FUNCT__
2268a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate"
2269d763cef2SBarry Smith /*@C
2270a6570f20SBarry Smith    TSMonitorLGCreate - Creates a line graph context for use with
2271d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
2272d763cef2SBarry Smith 
2273d763cef2SBarry Smith    Collective on TS
2274d763cef2SBarry Smith 
2275d763cef2SBarry Smith    Input Parameters:
2276d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
2277d763cef2SBarry Smith .  label - the title to put in the title bar
22787c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
2279d763cef2SBarry Smith -  m, n - the screen width and height in pixels
2280d763cef2SBarry Smith 
2281d763cef2SBarry Smith    Output Parameter:
2282d763cef2SBarry Smith .  draw - the drawing context
2283d763cef2SBarry Smith 
2284d763cef2SBarry Smith    Options Database Key:
2285a6570f20SBarry Smith .  -ts_monitor_draw - automatically sets line graph monitor
2286d763cef2SBarry Smith 
2287d763cef2SBarry Smith    Notes:
2288a6570f20SBarry Smith    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
2289d763cef2SBarry Smith 
2290d763cef2SBarry Smith    Level: intermediate
2291d763cef2SBarry Smith 
22927c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
2293d763cef2SBarry Smith 
2294a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet()
22957c922b88SBarry Smith 
2296d763cef2SBarry Smith @*/
2297a80ad3e0SBarry Smith PetscErrorCode  TSMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2298d763cef2SBarry Smith {
2299b0a32e0cSBarry Smith   PetscDraw      win;
2300dfbe8321SBarry Smith   PetscErrorCode ierr;
230170ddbef3SBarry Smith   PetscDrawAxis  axis;
2302d763cef2SBarry Smith 
2303d763cef2SBarry Smith   PetscFunctionBegin;
2304a80ad3e0SBarry Smith   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr);
2305b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
2306b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
2307b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
230870ddbef3SBarry Smith   ierr = PetscDrawLGGetAxis(*draw,&axis);CHKERRQ(ierr);
230970ddbef3SBarry Smith   ierr = PetscDrawAxisSetLabels(axis,"Integration Time","Iteration","Time");CHKERRQ(ierr);
231052e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
2311d763cef2SBarry Smith   PetscFunctionReturn(0);
2312d763cef2SBarry Smith }
2313d763cef2SBarry Smith 
23144a2ae208SSatish Balay #undef __FUNCT__
2315a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG"
2316a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
2317d763cef2SBarry Smith {
2318b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
231987828ca2SBarry Smith   PetscReal      x,y = ptime;
2320dfbe8321SBarry Smith   PetscErrorCode ierr;
2321d763cef2SBarry Smith 
2322d763cef2SBarry Smith   PetscFunctionBegin;
23237c922b88SBarry Smith   if (!monctx) {
23247c922b88SBarry Smith     MPI_Comm      comm;
2325b0a32e0cSBarry Smith     PetscViewer   viewer;
232670ddbef3SBarry Smith     PetscDrawAxis axis;
23277c922b88SBarry Smith 
23287c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
2329b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
2330b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
233170ddbef3SBarry Smith     ierr   = PetscDrawLGGetAxis(lg,&axis);CHKERRQ(ierr);
233270ddbef3SBarry Smith     ierr   = PetscDrawAxisSetLabels(axis,"Integration Time","Iteration","Time");CHKERRQ(ierr);
23337c922b88SBarry Smith   }
23347c922b88SBarry Smith 
2335b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
233687828ca2SBarry Smith   x = (PetscReal)n;
2337b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2338d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
2339b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2340d763cef2SBarry Smith   }
2341d763cef2SBarry Smith   PetscFunctionReturn(0);
2342d763cef2SBarry Smith }
2343d763cef2SBarry Smith 
23444a2ae208SSatish Balay #undef __FUNCT__
2345a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy"
2346d763cef2SBarry Smith /*@C
2347a6570f20SBarry Smith    TSMonitorLGDestroy - Destroys a line graph context that was created
2348a6570f20SBarry Smith    with TSMonitorLGCreate().
2349d763cef2SBarry Smith 
2350b0a32e0cSBarry Smith    Collective on PetscDrawLG
2351d763cef2SBarry Smith 
2352d763cef2SBarry Smith    Input Parameter:
2353d763cef2SBarry Smith .  draw - the drawing context
2354d763cef2SBarry Smith 
2355d763cef2SBarry Smith    Level: intermediate
2356d763cef2SBarry Smith 
2357d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
2358d763cef2SBarry Smith 
2359a6570f20SBarry Smith .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
2360d763cef2SBarry Smith @*/
23616bf464f9SBarry Smith PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
2362d763cef2SBarry Smith {
2363b0a32e0cSBarry Smith   PetscDraw      draw;
2364dfbe8321SBarry Smith   PetscErrorCode ierr;
2365d763cef2SBarry Smith 
2366d763cef2SBarry Smith   PetscFunctionBegin;
23676bf464f9SBarry Smith   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
23686bf464f9SBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
2369b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
2370d763cef2SBarry Smith   PetscFunctionReturn(0);
2371d763cef2SBarry Smith }
2372d763cef2SBarry Smith 
23734a2ae208SSatish Balay #undef __FUNCT__
23744a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
2375d763cef2SBarry Smith /*@
2376b8123daeSJed Brown    TSGetTime - Gets the time of the most recently completed step.
2377d763cef2SBarry Smith 
2378d763cef2SBarry Smith    Not Collective
2379d763cef2SBarry Smith 
2380d763cef2SBarry Smith    Input Parameter:
2381d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
2382d763cef2SBarry Smith 
2383d763cef2SBarry Smith    Output Parameter:
2384d763cef2SBarry Smith .  t  - the current time
2385d763cef2SBarry Smith 
2386d763cef2SBarry Smith    Level: beginner
2387d763cef2SBarry Smith 
2388b8123daeSJed Brown    Note:
2389b8123daeSJed Brown    When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
2390b8123daeSJed Brown    TSSetPreStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.
2391b8123daeSJed Brown 
2392d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
2393d763cef2SBarry Smith 
2394d763cef2SBarry Smith .keywords: TS, get, time
2395d763cef2SBarry Smith @*/
23967087cfbeSBarry Smith PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
2397d763cef2SBarry Smith {
2398d763cef2SBarry Smith   PetscFunctionBegin;
23990700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2400f7cf8827SBarry Smith   PetscValidRealPointer(t,2);
2401d763cef2SBarry Smith   *t = ts->ptime;
2402d763cef2SBarry Smith   PetscFunctionReturn(0);
2403d763cef2SBarry Smith }
2404d763cef2SBarry Smith 
24054a2ae208SSatish Balay #undef __FUNCT__
24066a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime"
24076a4d4014SLisandro Dalcin /*@
24086a4d4014SLisandro Dalcin    TSSetTime - Allows one to reset the time.
24096a4d4014SLisandro Dalcin 
24103f9fe445SBarry Smith    Logically Collective on TS
24116a4d4014SLisandro Dalcin 
24126a4d4014SLisandro Dalcin    Input Parameters:
24136a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
24146a4d4014SLisandro Dalcin -  time - the time
24156a4d4014SLisandro Dalcin 
24166a4d4014SLisandro Dalcin    Level: intermediate
24176a4d4014SLisandro Dalcin 
24186a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration()
24196a4d4014SLisandro Dalcin 
24206a4d4014SLisandro Dalcin .keywords: TS, set, time
24216a4d4014SLisandro Dalcin @*/
24227087cfbeSBarry Smith PetscErrorCode  TSSetTime(TS ts, PetscReal t)
24236a4d4014SLisandro Dalcin {
24246a4d4014SLisandro Dalcin   PetscFunctionBegin;
24250700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2426c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,t,2);
24276a4d4014SLisandro Dalcin   ts->ptime = t;
24286a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
24296a4d4014SLisandro Dalcin }
24306a4d4014SLisandro Dalcin 
24316a4d4014SLisandro Dalcin #undef __FUNCT__
24324a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
2433d763cef2SBarry Smith /*@C
2434d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
2435d763cef2SBarry Smith    TS options in the database.
2436d763cef2SBarry Smith 
24373f9fe445SBarry Smith    Logically Collective on TS
2438d763cef2SBarry Smith 
2439d763cef2SBarry Smith    Input Parameter:
2440d763cef2SBarry Smith +  ts     - The TS context
2441d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
2442d763cef2SBarry Smith 
2443d763cef2SBarry Smith    Notes:
2444d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
2445d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
2446d763cef2SBarry Smith    hyphen.
2447d763cef2SBarry Smith 
2448d763cef2SBarry Smith    Level: advanced
2449d763cef2SBarry Smith 
2450d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
2451d763cef2SBarry Smith 
2452d763cef2SBarry Smith .seealso: TSSetFromOptions()
2453d763cef2SBarry Smith 
2454d763cef2SBarry Smith @*/
24557087cfbeSBarry Smith PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2456d763cef2SBarry Smith {
2457dfbe8321SBarry Smith   PetscErrorCode ierr;
2458089b2837SJed Brown   SNES           snes;
2459d763cef2SBarry Smith 
2460d763cef2SBarry Smith   PetscFunctionBegin;
24610700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2462d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2463089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2464089b2837SJed Brown   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2465d763cef2SBarry Smith   PetscFunctionReturn(0);
2466d763cef2SBarry Smith }
2467d763cef2SBarry Smith 
2468d763cef2SBarry Smith 
24694a2ae208SSatish Balay #undef __FUNCT__
24704a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
2471d763cef2SBarry Smith /*@C
2472d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2473d763cef2SBarry Smith    TS options in the database.
2474d763cef2SBarry Smith 
24753f9fe445SBarry Smith    Logically Collective on TS
2476d763cef2SBarry Smith 
2477d763cef2SBarry Smith    Input Parameter:
2478d763cef2SBarry Smith +  ts     - The TS context
2479d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
2480d763cef2SBarry Smith 
2481d763cef2SBarry Smith    Notes:
2482d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
2483d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
2484d763cef2SBarry Smith    hyphen.
2485d763cef2SBarry Smith 
2486d763cef2SBarry Smith    Level: advanced
2487d763cef2SBarry Smith 
2488d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
2489d763cef2SBarry Smith 
2490d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
2491d763cef2SBarry Smith 
2492d763cef2SBarry Smith @*/
24937087cfbeSBarry Smith PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2494d763cef2SBarry Smith {
2495dfbe8321SBarry Smith   PetscErrorCode ierr;
2496089b2837SJed Brown   SNES           snes;
2497d763cef2SBarry Smith 
2498d763cef2SBarry Smith   PetscFunctionBegin;
24990700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2500d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2501089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2502089b2837SJed Brown   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2503d763cef2SBarry Smith   PetscFunctionReturn(0);
2504d763cef2SBarry Smith }
2505d763cef2SBarry Smith 
25064a2ae208SSatish Balay #undef __FUNCT__
25074a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
2508d763cef2SBarry Smith /*@C
2509d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
2510d763cef2SBarry Smith    TS options in the database.
2511d763cef2SBarry Smith 
2512d763cef2SBarry Smith    Not Collective
2513d763cef2SBarry Smith 
2514d763cef2SBarry Smith    Input Parameter:
2515d763cef2SBarry Smith .  ts - The TS context
2516d763cef2SBarry Smith 
2517d763cef2SBarry Smith    Output Parameter:
2518d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
2519d763cef2SBarry Smith 
2520d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
2521d763cef2SBarry Smith    sufficient length to hold the prefix.
2522d763cef2SBarry Smith 
2523d763cef2SBarry Smith    Level: intermediate
2524d763cef2SBarry Smith 
2525d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
2526d763cef2SBarry Smith 
2527d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
2528d763cef2SBarry Smith @*/
25297087cfbeSBarry Smith PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2530d763cef2SBarry Smith {
2531dfbe8321SBarry Smith   PetscErrorCode ierr;
2532d763cef2SBarry Smith 
2533d763cef2SBarry Smith   PetscFunctionBegin;
25340700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
25354482741eSBarry Smith   PetscValidPointer(prefix,2);
2536d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2537d763cef2SBarry Smith   PetscFunctionReturn(0);
2538d763cef2SBarry Smith }
2539d763cef2SBarry Smith 
25404a2ae208SSatish Balay #undef __FUNCT__
25414a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
2542d763cef2SBarry Smith /*@C
2543d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2544d763cef2SBarry Smith 
2545d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
2546d763cef2SBarry Smith 
2547d763cef2SBarry Smith    Input Parameter:
2548d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
2549d763cef2SBarry Smith 
2550d763cef2SBarry Smith    Output Parameters:
2551d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
2552d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
2553089b2837SJed Brown .  func - Function to compute the Jacobian of the RHS
2554d763cef2SBarry Smith -  ctx - User-defined context for Jacobian evaluation routine
2555d763cef2SBarry Smith 
2556d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2557d763cef2SBarry Smith 
2558d763cef2SBarry Smith    Level: intermediate
2559d763cef2SBarry Smith 
256026d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2561d763cef2SBarry Smith 
2562d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
2563d763cef2SBarry Smith @*/
2564089b2837SJed Brown PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2565d763cef2SBarry Smith {
2566089b2837SJed Brown   PetscErrorCode ierr;
2567089b2837SJed Brown   SNES           snes;
256824989b8cSPeter Brune   DM             dm;
2569089b2837SJed Brown 
2570d763cef2SBarry Smith   PetscFunctionBegin;
2571089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2572089b2837SJed Brown   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
257324989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
257424989b8cSPeter Brune   ierr = DMTSGetRHSJacobian(dm,func,ctx);CHKERRQ(ierr);
2575d763cef2SBarry Smith   PetscFunctionReturn(0);
2576d763cef2SBarry Smith }
2577d763cef2SBarry Smith 
25781713a123SBarry Smith #undef __FUNCT__
25792eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian"
25802eca1d9cSJed Brown /*@C
25812eca1d9cSJed Brown    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
25822eca1d9cSJed Brown 
25832eca1d9cSJed Brown    Not Collective, but parallel objects are returned if TS is parallel
25842eca1d9cSJed Brown 
25852eca1d9cSJed Brown    Input Parameter:
25862eca1d9cSJed Brown .  ts  - The TS context obtained from TSCreate()
25872eca1d9cSJed Brown 
25882eca1d9cSJed Brown    Output Parameters:
25892eca1d9cSJed Brown +  A   - The Jacobian of F(t,U,U_t)
25902eca1d9cSJed Brown .  B   - The preconditioner matrix, often the same as A
25912eca1d9cSJed Brown .  f   - The function to compute the matrices
25922eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine
25932eca1d9cSJed Brown 
25942eca1d9cSJed Brown    Notes: You can pass in PETSC_NULL for any return argument you do not need.
25952eca1d9cSJed Brown 
25962eca1d9cSJed Brown    Level: advanced
25972eca1d9cSJed Brown 
25982eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
25992eca1d9cSJed Brown 
26002eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian
26012eca1d9cSJed Brown @*/
26027087cfbeSBarry Smith PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
26032eca1d9cSJed Brown {
2604089b2837SJed Brown   PetscErrorCode ierr;
2605089b2837SJed Brown   SNES           snes;
260624989b8cSPeter Brune   DM             dm;
26072eca1d9cSJed Brown   PetscFunctionBegin;
2608089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2609089b2837SJed Brown   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
261024989b8cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
261124989b8cSPeter Brune   ierr = DMTSGetIJacobian(dm,f,ctx);CHKERRQ(ierr);
26122eca1d9cSJed Brown   PetscFunctionReturn(0);
26132eca1d9cSJed Brown }
26142eca1d9cSJed Brown 
26156083293cSBarry Smith typedef struct {
26166083293cSBarry Smith   PetscViewer viewer;
26176083293cSBarry Smith   Vec         initialsolution;
26186083293cSBarry Smith   PetscBool   showinitial;
26196083293cSBarry Smith } TSMonitorSolutionCtx;
26206083293cSBarry Smith 
26212eca1d9cSJed Brown #undef __FUNCT__
2622a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution"
26231713a123SBarry Smith /*@C
2624a6570f20SBarry Smith    TSMonitorSolution - Monitors progress of the TS solvers by calling
26251713a123SBarry Smith    VecView() for the solution at each timestep
26261713a123SBarry Smith 
26271713a123SBarry Smith    Collective on TS
26281713a123SBarry Smith 
26291713a123SBarry Smith    Input Parameters:
26301713a123SBarry Smith +  ts - the TS context
26311713a123SBarry Smith .  step - current time-step
2632142b95e3SSatish Balay .  ptime - current time
26331713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
26341713a123SBarry Smith 
26351713a123SBarry Smith    Level: intermediate
26361713a123SBarry Smith 
26371713a123SBarry Smith .keywords: TS,  vector, monitor, view
26381713a123SBarry Smith 
2639a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
26401713a123SBarry Smith @*/
26417087cfbeSBarry Smith PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
26421713a123SBarry Smith {
2643dfbe8321SBarry Smith   PetscErrorCode       ierr;
26446083293cSBarry Smith   TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy;
26451713a123SBarry Smith 
26461713a123SBarry Smith   PetscFunctionBegin;
26476083293cSBarry Smith   if (!step && ictx->showinitial) {
26486083293cSBarry Smith     if (!ictx->initialsolution) {
26496083293cSBarry Smith       ierr = VecDuplicate(x,&ictx->initialsolution);CHKERRQ(ierr);
26501713a123SBarry Smith     }
26516083293cSBarry Smith     ierr = VecCopy(x,ictx->initialsolution);CHKERRQ(ierr);
26526083293cSBarry Smith   }
26536083293cSBarry Smith   if (ictx->showinitial) {
26546083293cSBarry Smith     PetscReal pause;
26556083293cSBarry Smith     ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr);
26566083293cSBarry Smith     ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr);
26576083293cSBarry Smith     ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr);
26586083293cSBarry Smith     ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr);
26596083293cSBarry Smith     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr);
26606083293cSBarry Smith   }
26616083293cSBarry Smith   ierr = VecView(x,ictx->viewer);CHKERRQ(ierr);
26626083293cSBarry Smith   if (ictx->showinitial) {
26636083293cSBarry Smith     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr);
26646083293cSBarry Smith   }
26651713a123SBarry Smith   PetscFunctionReturn(0);
26661713a123SBarry Smith }
26671713a123SBarry Smith 
26681713a123SBarry Smith 
26696c699258SBarry Smith #undef __FUNCT__
26706083293cSBarry Smith #define __FUNCT__ "TSMonitorSolutionDestroy"
26716083293cSBarry Smith /*@C
26726083293cSBarry Smith    TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution
26736083293cSBarry Smith 
26746083293cSBarry Smith    Collective on TS
26756083293cSBarry Smith 
26766083293cSBarry Smith    Input Parameters:
26776083293cSBarry Smith .    ctx - the monitor context
26786083293cSBarry Smith 
26796083293cSBarry Smith    Level: intermediate
26806083293cSBarry Smith 
26816083293cSBarry Smith .keywords: TS,  vector, monitor, view
26826083293cSBarry Smith 
26836083293cSBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
26846083293cSBarry Smith @*/
26856083293cSBarry Smith PetscErrorCode  TSMonitorSolutionDestroy(void **ctx)
26866083293cSBarry Smith {
26876083293cSBarry Smith   PetscErrorCode       ierr;
26886083293cSBarry Smith   TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx;
26896083293cSBarry Smith 
26906083293cSBarry Smith   PetscFunctionBegin;
26916083293cSBarry Smith   ierr = PetscViewerDestroy(&ictx->viewer);CHKERRQ(ierr);
26926083293cSBarry Smith   ierr = VecDestroy(&ictx->initialsolution);CHKERRQ(ierr);
26936083293cSBarry Smith   ierr = PetscFree(ictx);CHKERRQ(ierr);
26946083293cSBarry Smith   PetscFunctionReturn(0);
26956083293cSBarry Smith }
26966083293cSBarry Smith 
26976083293cSBarry Smith #undef __FUNCT__
26986083293cSBarry Smith #define __FUNCT__ "TSMonitorSolutionCreate"
26996083293cSBarry Smith /*@C
27006083293cSBarry Smith    TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution
27016083293cSBarry Smith 
27026083293cSBarry Smith    Collective on TS
27036083293cSBarry Smith 
27046083293cSBarry Smith    Input Parameter:
27056083293cSBarry Smith .    ts - time-step context
27066083293cSBarry Smith 
27076083293cSBarry Smith    Output Patameter:
27086083293cSBarry Smith .    ctx - the monitor context
27096083293cSBarry Smith 
27106083293cSBarry Smith    Level: intermediate
27116083293cSBarry Smith 
27126083293cSBarry Smith .keywords: TS,  vector, monitor, view
27136083293cSBarry Smith 
27146083293cSBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
27156083293cSBarry Smith @*/
27166083293cSBarry Smith PetscErrorCode  TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx)
27176083293cSBarry Smith {
27186083293cSBarry Smith   PetscErrorCode       ierr;
27196083293cSBarry Smith   TSMonitorSolutionCtx *ictx;
27206083293cSBarry Smith 
27216083293cSBarry Smith   PetscFunctionBegin;
27226083293cSBarry Smith   ierr = PetscNew(TSMonitorSolutionCtx,&ictx);CHKERRQ(ierr);
27236083293cSBarry Smith   *ctx = (void*)ictx;
27246083293cSBarry Smith   if (!viewer) {
27256083293cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
27266083293cSBarry Smith   }
27276083293cSBarry Smith   ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
27286083293cSBarry Smith   ictx->viewer      = viewer;
27296083293cSBarry Smith   ictx->showinitial = PETSC_FALSE;
27306083293cSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);CHKERRQ(ierr);
27316083293cSBarry Smith   PetscFunctionReturn(0);
27326083293cSBarry Smith }
27336083293cSBarry Smith 
27346083293cSBarry Smith #undef __FUNCT__
27356c699258SBarry Smith #define __FUNCT__ "TSSetDM"
27366c699258SBarry Smith /*@
27376c699258SBarry Smith    TSSetDM - Sets the DM that may be used by some preconditioners
27386c699258SBarry Smith 
27393f9fe445SBarry Smith    Logically Collective on TS and DM
27406c699258SBarry Smith 
27416c699258SBarry Smith    Input Parameters:
27426c699258SBarry Smith +  ts - the preconditioner context
27436c699258SBarry Smith -  dm - the dm
27446c699258SBarry Smith 
27456c699258SBarry Smith    Level: intermediate
27466c699258SBarry Smith 
27476c699258SBarry Smith 
27486c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
27496c699258SBarry Smith @*/
27507087cfbeSBarry Smith PetscErrorCode  TSSetDM(TS ts,DM dm)
27516c699258SBarry Smith {
27526c699258SBarry Smith   PetscErrorCode ierr;
2753089b2837SJed Brown   SNES           snes;
275424989b8cSPeter Brune   TSDM           tsdm;
27556c699258SBarry Smith 
27566c699258SBarry Smith   PetscFunctionBegin;
27570700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
275870663e4aSLisandro Dalcin   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
275924989b8cSPeter Brune   if (ts->dm) {               /* Move the TSDM context over to the new DM unless the new DM already has one */
276024989b8cSPeter Brune     PetscContainer oldcontainer,container;
276124989b8cSPeter Brune     ierr = PetscObjectQuery((PetscObject)ts->dm,"TSDM",(PetscObject*)&oldcontainer);CHKERRQ(ierr);
276224989b8cSPeter Brune     ierr = PetscObjectQuery((PetscObject)dm,"TSDM",(PetscObject*)&container);CHKERRQ(ierr);
276324989b8cSPeter Brune     if (oldcontainer && !container) {
276424989b8cSPeter Brune       ierr = DMTSCopyContext(ts->dm,dm);CHKERRQ(ierr);
276524989b8cSPeter Brune       ierr = DMTSGetContext(ts->dm,&tsdm);CHKERRQ(ierr);
276624989b8cSPeter Brune       if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */
276724989b8cSPeter Brune         tsdm->originaldm = dm;
276824989b8cSPeter Brune       }
276924989b8cSPeter Brune     }
27706bf464f9SBarry Smith     ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
277124989b8cSPeter Brune   }
27726c699258SBarry Smith   ts->dm = dm;
2773089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2774089b2837SJed Brown   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
27756c699258SBarry Smith   PetscFunctionReturn(0);
27766c699258SBarry Smith }
27776c699258SBarry Smith 
27786c699258SBarry Smith #undef __FUNCT__
27796c699258SBarry Smith #define __FUNCT__ "TSGetDM"
27806c699258SBarry Smith /*@
27816c699258SBarry Smith    TSGetDM - Gets the DM that may be used by some preconditioners
27826c699258SBarry Smith 
27833f9fe445SBarry Smith    Not Collective
27846c699258SBarry Smith 
27856c699258SBarry Smith    Input Parameter:
27866c699258SBarry Smith . ts - the preconditioner context
27876c699258SBarry Smith 
27886c699258SBarry Smith    Output Parameter:
27896c699258SBarry Smith .  dm - the dm
27906c699258SBarry Smith 
27916c699258SBarry Smith    Level: intermediate
27926c699258SBarry Smith 
27936c699258SBarry Smith 
27946c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
27956c699258SBarry Smith @*/
27967087cfbeSBarry Smith PetscErrorCode  TSGetDM(TS ts,DM *dm)
27976c699258SBarry Smith {
2798496e6a7aSJed Brown   PetscErrorCode ierr;
2799496e6a7aSJed Brown 
28006c699258SBarry Smith   PetscFunctionBegin;
28010700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2802496e6a7aSJed Brown   if (!ts->dm) {
2803496e6a7aSJed Brown     ierr = DMShellCreate(((PetscObject)ts)->comm,&ts->dm);CHKERRQ(ierr);
2804496e6a7aSJed Brown     if (ts->snes) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
2805496e6a7aSJed Brown   }
28066c699258SBarry Smith   *dm = ts->dm;
28076c699258SBarry Smith   PetscFunctionReturn(0);
28086c699258SBarry Smith }
28091713a123SBarry Smith 
28100f5c6efeSJed Brown #undef __FUNCT__
28110f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction"
28120f5c6efeSJed Brown /*@
28130f5c6efeSJed Brown    SNESTSFormFunction - Function to evaluate nonlinear residual
28140f5c6efeSJed Brown 
28153f9fe445SBarry Smith    Logically Collective on SNES
28160f5c6efeSJed Brown 
28170f5c6efeSJed Brown    Input Parameter:
2818d42a1c89SJed Brown + snes - nonlinear solver
28190f5c6efeSJed Brown . X - the current state at which to evaluate the residual
2820d42a1c89SJed Brown - ctx - user context, must be a TS
28210f5c6efeSJed Brown 
28220f5c6efeSJed Brown    Output Parameter:
28230f5c6efeSJed Brown . F - the nonlinear residual
28240f5c6efeSJed Brown 
28250f5c6efeSJed Brown    Notes:
28260f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
28270f5c6efeSJed Brown    It is most frequently passed to MatFDColoringSetFunction().
28280f5c6efeSJed Brown 
28290f5c6efeSJed Brown    Level: advanced
28300f5c6efeSJed Brown 
28310f5c6efeSJed Brown .seealso: SNESSetFunction(), MatFDColoringSetFunction()
28320f5c6efeSJed Brown @*/
28337087cfbeSBarry Smith PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
28340f5c6efeSJed Brown {
28350f5c6efeSJed Brown   TS ts = (TS)ctx;
28360f5c6efeSJed Brown   PetscErrorCode ierr;
28370f5c6efeSJed Brown 
28380f5c6efeSJed Brown   PetscFunctionBegin;
28390f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
28400f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
28410f5c6efeSJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
28420f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
28430f5c6efeSJed Brown   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
28440f5c6efeSJed Brown   PetscFunctionReturn(0);
28450f5c6efeSJed Brown }
28460f5c6efeSJed Brown 
28470f5c6efeSJed Brown #undef __FUNCT__
28480f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian"
28490f5c6efeSJed Brown /*@
28500f5c6efeSJed Brown    SNESTSFormJacobian - Function to evaluate the Jacobian
28510f5c6efeSJed Brown 
28520f5c6efeSJed Brown    Collective on SNES
28530f5c6efeSJed Brown 
28540f5c6efeSJed Brown    Input Parameter:
28550f5c6efeSJed Brown + snes - nonlinear solver
28560f5c6efeSJed Brown . X - the current state at which to evaluate the residual
28570f5c6efeSJed Brown - ctx - user context, must be a TS
28580f5c6efeSJed Brown 
28590f5c6efeSJed Brown    Output Parameter:
28600f5c6efeSJed Brown + A - the Jacobian
28610f5c6efeSJed Brown . B - the preconditioning matrix (may be the same as A)
28620f5c6efeSJed Brown - flag - indicates any structure change in the matrix
28630f5c6efeSJed Brown 
28640f5c6efeSJed Brown    Notes:
28650f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
28660f5c6efeSJed Brown 
28670f5c6efeSJed Brown    Level: developer
28680f5c6efeSJed Brown 
28690f5c6efeSJed Brown .seealso: SNESSetJacobian()
28700f5c6efeSJed Brown @*/
28717087cfbeSBarry Smith PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
28720f5c6efeSJed Brown {
28730f5c6efeSJed Brown   TS ts = (TS)ctx;
28740f5c6efeSJed Brown   PetscErrorCode ierr;
28750f5c6efeSJed Brown 
28760f5c6efeSJed Brown   PetscFunctionBegin;
28770f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
28780f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
28790f5c6efeSJed Brown   PetscValidPointer(A,3);
28800f5c6efeSJed Brown   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
28810f5c6efeSJed Brown   PetscValidPointer(B,4);
28820f5c6efeSJed Brown   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
28830f5c6efeSJed Brown   PetscValidPointer(flag,5);
28840f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
28850f5c6efeSJed Brown   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
28860f5c6efeSJed Brown   PetscFunctionReturn(0);
28870f5c6efeSJed Brown }
2888325fc9f4SBarry Smith 
28890e4ef248SJed Brown #undef __FUNCT__
28900e4ef248SJed Brown #define __FUNCT__ "TSComputeRHSFunctionLinear"
28910e4ef248SJed Brown /*@C
28920e4ef248SJed Brown    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only
28930e4ef248SJed Brown 
28940e4ef248SJed Brown    Collective on TS
28950e4ef248SJed Brown 
28960e4ef248SJed Brown    Input Arguments:
28970e4ef248SJed Brown +  ts - time stepping context
28980e4ef248SJed Brown .  t - time at which to evaluate
28990e4ef248SJed Brown .  X - state at which to evaluate
29000e4ef248SJed Brown -  ctx - context
29010e4ef248SJed Brown 
29020e4ef248SJed Brown    Output Arguments:
29030e4ef248SJed Brown .  F - right hand side
29040e4ef248SJed Brown 
29050e4ef248SJed Brown    Level: intermediate
29060e4ef248SJed Brown 
29070e4ef248SJed Brown    Notes:
29080e4ef248SJed Brown    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
29090e4ef248SJed Brown    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
29100e4ef248SJed Brown 
29110e4ef248SJed Brown .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
29120e4ef248SJed Brown @*/
29130e4ef248SJed Brown PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
29140e4ef248SJed Brown {
29150e4ef248SJed Brown   PetscErrorCode ierr;
29160e4ef248SJed Brown   Mat Arhs,Brhs;
29170e4ef248SJed Brown   MatStructure flg2;
29180e4ef248SJed Brown 
29190e4ef248SJed Brown   PetscFunctionBegin;
29200e4ef248SJed Brown   ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
29210e4ef248SJed Brown   ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
29220e4ef248SJed Brown   ierr = MatMult(Arhs,X,F);CHKERRQ(ierr);
29230e4ef248SJed Brown   PetscFunctionReturn(0);
29240e4ef248SJed Brown }
29250e4ef248SJed Brown 
29260e4ef248SJed Brown #undef __FUNCT__
29270e4ef248SJed Brown #define __FUNCT__ "TSComputeRHSJacobianConstant"
29280e4ef248SJed Brown /*@C
29290e4ef248SJed Brown    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
29300e4ef248SJed Brown 
29310e4ef248SJed Brown    Collective on TS
29320e4ef248SJed Brown 
29330e4ef248SJed Brown    Input Arguments:
29340e4ef248SJed Brown +  ts - time stepping context
29350e4ef248SJed Brown .  t - time at which to evaluate
29360e4ef248SJed Brown .  X - state at which to evaluate
29370e4ef248SJed Brown -  ctx - context
29380e4ef248SJed Brown 
29390e4ef248SJed Brown    Output Arguments:
29400e4ef248SJed Brown +  A - pointer to operator
29410e4ef248SJed Brown .  B - pointer to preconditioning matrix
29420e4ef248SJed Brown -  flg - matrix structure flag
29430e4ef248SJed Brown 
29440e4ef248SJed Brown    Level: intermediate
29450e4ef248SJed Brown 
29460e4ef248SJed Brown    Notes:
29470e4ef248SJed Brown    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
29480e4ef248SJed Brown 
29490e4ef248SJed Brown .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
29500e4ef248SJed Brown @*/
29510e4ef248SJed Brown PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
29520e4ef248SJed Brown {
29530e4ef248SJed Brown 
29540e4ef248SJed Brown   PetscFunctionBegin;
29550e4ef248SJed Brown   *flg = SAME_PRECONDITIONER;
29560e4ef248SJed Brown   PetscFunctionReturn(0);
29570e4ef248SJed Brown }
29580e4ef248SJed Brown 
29590026cea9SSean Farley #undef __FUNCT__
29600026cea9SSean Farley #define __FUNCT__ "TSComputeIFunctionLinear"
29610026cea9SSean Farley /*@C
29620026cea9SSean Farley    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
29630026cea9SSean Farley 
29640026cea9SSean Farley    Collective on TS
29650026cea9SSean Farley 
29660026cea9SSean Farley    Input Arguments:
29670026cea9SSean Farley +  ts - time stepping context
29680026cea9SSean Farley .  t - time at which to evaluate
29690026cea9SSean Farley .  X - state at which to evaluate
29700026cea9SSean Farley .  Xdot - time derivative of state vector
29710026cea9SSean Farley -  ctx - context
29720026cea9SSean Farley 
29730026cea9SSean Farley    Output Arguments:
29740026cea9SSean Farley .  F - left hand side
29750026cea9SSean Farley 
29760026cea9SSean Farley    Level: intermediate
29770026cea9SSean Farley 
29780026cea9SSean Farley    Notes:
29790026cea9SSean 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
29800026cea9SSean Farley    user is required to write their own TSComputeIFunction.
29810026cea9SSean Farley    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
29820026cea9SSean Farley    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
29830026cea9SSean Farley 
29840026cea9SSean Farley .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
29850026cea9SSean Farley @*/
29860026cea9SSean Farley PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
29870026cea9SSean Farley {
29880026cea9SSean Farley   PetscErrorCode ierr;
29890026cea9SSean Farley   Mat A,B;
29900026cea9SSean Farley   MatStructure flg2;
29910026cea9SSean Farley 
29920026cea9SSean Farley   PetscFunctionBegin;
29930026cea9SSean Farley   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
29940026cea9SSean Farley   ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr);
29950026cea9SSean Farley   ierr = MatMult(A,Xdot,F);CHKERRQ(ierr);
29960026cea9SSean Farley   PetscFunctionReturn(0);
29970026cea9SSean Farley }
29980026cea9SSean Farley 
29990026cea9SSean Farley #undef __FUNCT__
30000026cea9SSean Farley #define __FUNCT__ "TSComputeIJacobianConstant"
30010026cea9SSean Farley /*@C
300224e94211SJed Brown    TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent.
30030026cea9SSean Farley 
30040026cea9SSean Farley    Collective on TS
30050026cea9SSean Farley 
30060026cea9SSean Farley    Input Arguments:
30070026cea9SSean Farley +  ts - time stepping context
30080026cea9SSean Farley .  t - time at which to evaluate
30090026cea9SSean Farley .  X - state at which to evaluate
30100026cea9SSean Farley .  Xdot - time derivative of state vector
30110026cea9SSean Farley .  shift - shift to apply
30120026cea9SSean Farley -  ctx - context
30130026cea9SSean Farley 
30140026cea9SSean Farley    Output Arguments:
30150026cea9SSean Farley +  A - pointer to operator
30160026cea9SSean Farley .  B - pointer to preconditioning matrix
30170026cea9SSean Farley -  flg - matrix structure flag
30180026cea9SSean Farley 
30190026cea9SSean Farley    Level: intermediate
30200026cea9SSean Farley 
30210026cea9SSean Farley    Notes:
30220026cea9SSean Farley    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
30230026cea9SSean Farley 
30240026cea9SSean Farley .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
30250026cea9SSean Farley @*/
30260026cea9SSean Farley PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
30270026cea9SSean Farley {
30280026cea9SSean Farley 
30290026cea9SSean Farley   PetscFunctionBegin;
30300026cea9SSean Farley   *flg = SAME_PRECONDITIONER;
30310026cea9SSean Farley   PetscFunctionReturn(0);
30320026cea9SSean Farley }
30330026cea9SSean Farley 
30344af1b03aSJed Brown 
30354af1b03aSJed Brown #undef __FUNCT__
30364af1b03aSJed Brown #define __FUNCT__ "TSGetConvergedReason"
30374af1b03aSJed Brown /*@
30384af1b03aSJed Brown    TSGetConvergedReason - Gets the reason the TS iteration was stopped.
30394af1b03aSJed Brown 
30404af1b03aSJed Brown    Not Collective
30414af1b03aSJed Brown 
30424af1b03aSJed Brown    Input Parameter:
30434af1b03aSJed Brown .  ts - the TS context
30444af1b03aSJed Brown 
30454af1b03aSJed Brown    Output Parameter:
30464af1b03aSJed Brown .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
30474af1b03aSJed Brown             manual pages for the individual convergence tests for complete lists
30484af1b03aSJed Brown 
30494af1b03aSJed Brown    Level: intermediate
30504af1b03aSJed Brown 
3051cd652676SJed Brown    Notes:
3052cd652676SJed Brown    Can only be called after the call to TSSolve() is complete.
30534af1b03aSJed Brown 
30544af1b03aSJed Brown .keywords: TS, nonlinear, set, convergence, test
30554af1b03aSJed Brown 
30564af1b03aSJed Brown .seealso: TSSetConvergenceTest(), TSConvergedReason
30574af1b03aSJed Brown @*/
30584af1b03aSJed Brown PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
30594af1b03aSJed Brown {
30604af1b03aSJed Brown   PetscFunctionBegin;
30614af1b03aSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
30624af1b03aSJed Brown   PetscValidPointer(reason,2);
30634af1b03aSJed Brown   *reason = ts->reason;
30644af1b03aSJed Brown   PetscFunctionReturn(0);
30654af1b03aSJed Brown }
30664af1b03aSJed Brown 
3067fb1732b5SBarry Smith #undef __FUNCT__
30685ef26d82SJed Brown #define __FUNCT__ "TSGetSNESIterations"
30699f67acb7SJed Brown /*@
30705ef26d82SJed Brown    TSGetSNESIterations - Gets the total number of nonlinear iterations
30719f67acb7SJed Brown    used by the time integrator.
30729f67acb7SJed Brown 
30739f67acb7SJed Brown    Not Collective
30749f67acb7SJed Brown 
30759f67acb7SJed Brown    Input Parameter:
30769f67acb7SJed Brown .  ts - TS context
30779f67acb7SJed Brown 
30789f67acb7SJed Brown    Output Parameter:
30799f67acb7SJed Brown .  nits - number of nonlinear iterations
30809f67acb7SJed Brown 
30819f67acb7SJed Brown    Notes:
30829f67acb7SJed Brown    This counter is reset to zero for each successive call to TSSolve().
30839f67acb7SJed Brown 
30849f67acb7SJed Brown    Level: intermediate
30859f67acb7SJed Brown 
30869f67acb7SJed Brown .keywords: TS, get, number, nonlinear, iterations
30879f67acb7SJed Brown 
30885ef26d82SJed Brown .seealso:  TSGetKSPIterations()
30899f67acb7SJed Brown @*/
30905ef26d82SJed Brown PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
30919f67acb7SJed Brown {
30929f67acb7SJed Brown   PetscFunctionBegin;
30939f67acb7SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
30949f67acb7SJed Brown   PetscValidIntPointer(nits,2);
30955ef26d82SJed Brown   *nits = ts->snes_its;
30969f67acb7SJed Brown   PetscFunctionReturn(0);
30979f67acb7SJed Brown }
30989f67acb7SJed Brown 
30999f67acb7SJed Brown #undef __FUNCT__
31005ef26d82SJed Brown #define __FUNCT__ "TSGetKSPIterations"
31019f67acb7SJed Brown /*@
31025ef26d82SJed Brown    TSGetKSPIterations - Gets the total number of linear iterations
31039f67acb7SJed Brown    used by the time integrator.
31049f67acb7SJed Brown 
31059f67acb7SJed Brown    Not Collective
31069f67acb7SJed Brown 
31079f67acb7SJed Brown    Input Parameter:
31089f67acb7SJed Brown .  ts - TS context
31099f67acb7SJed Brown 
31109f67acb7SJed Brown    Output Parameter:
31119f67acb7SJed Brown .  lits - number of linear iterations
31129f67acb7SJed Brown 
31139f67acb7SJed Brown    Notes:
31149f67acb7SJed Brown    This counter is reset to zero for each successive call to TSSolve().
31159f67acb7SJed Brown 
31169f67acb7SJed Brown    Level: intermediate
31179f67acb7SJed Brown 
31189f67acb7SJed Brown .keywords: TS, get, number, linear, iterations
31199f67acb7SJed Brown 
31205ef26d82SJed Brown .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
31219f67acb7SJed Brown @*/
31225ef26d82SJed Brown PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
31239f67acb7SJed Brown {
31249f67acb7SJed Brown   PetscFunctionBegin;
31259f67acb7SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
31269f67acb7SJed Brown   PetscValidIntPointer(lits,2);
31275ef26d82SJed Brown   *lits = ts->ksp_its;
31289f67acb7SJed Brown   PetscFunctionReturn(0);
31299f67acb7SJed Brown }
31309f67acb7SJed Brown 
31319f67acb7SJed Brown #undef __FUNCT__
3132cef5090cSJed Brown #define __FUNCT__ "TSGetStepRejections"
3133cef5090cSJed Brown /*@
3134cef5090cSJed Brown    TSGetStepRejections - Gets the total number of rejected steps.
3135cef5090cSJed Brown 
3136cef5090cSJed Brown    Not Collective
3137cef5090cSJed Brown 
3138cef5090cSJed Brown    Input Parameter:
3139cef5090cSJed Brown .  ts - TS context
3140cef5090cSJed Brown 
3141cef5090cSJed Brown    Output Parameter:
3142cef5090cSJed Brown .  rejects - number of steps rejected
3143cef5090cSJed Brown 
3144cef5090cSJed Brown    Notes:
3145cef5090cSJed Brown    This counter is reset to zero for each successive call to TSSolve().
3146cef5090cSJed Brown 
3147cef5090cSJed Brown    Level: intermediate
3148cef5090cSJed Brown 
3149cef5090cSJed Brown .keywords: TS, get, number
3150cef5090cSJed Brown 
31515ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
3152cef5090cSJed Brown @*/
3153cef5090cSJed Brown PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
3154cef5090cSJed Brown {
3155cef5090cSJed Brown   PetscFunctionBegin;
3156cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3157cef5090cSJed Brown   PetscValidIntPointer(rejects,2);
3158cef5090cSJed Brown   *rejects = ts->reject;
3159cef5090cSJed Brown   PetscFunctionReturn(0);
3160cef5090cSJed Brown }
3161cef5090cSJed Brown 
3162cef5090cSJed Brown #undef __FUNCT__
3163cef5090cSJed Brown #define __FUNCT__ "TSGetSNESFailures"
3164cef5090cSJed Brown /*@
3165cef5090cSJed Brown    TSGetSNESFailures - Gets the total number of failed SNES solves
3166cef5090cSJed Brown 
3167cef5090cSJed Brown    Not Collective
3168cef5090cSJed Brown 
3169cef5090cSJed Brown    Input Parameter:
3170cef5090cSJed Brown .  ts - TS context
3171cef5090cSJed Brown 
3172cef5090cSJed Brown    Output Parameter:
3173cef5090cSJed Brown .  fails - number of failed nonlinear solves
3174cef5090cSJed Brown 
3175cef5090cSJed Brown    Notes:
3176cef5090cSJed Brown    This counter is reset to zero for each successive call to TSSolve().
3177cef5090cSJed Brown 
3178cef5090cSJed Brown    Level: intermediate
3179cef5090cSJed Brown 
3180cef5090cSJed Brown .keywords: TS, get, number
3181cef5090cSJed Brown 
31825ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
3183cef5090cSJed Brown @*/
3184cef5090cSJed Brown PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
3185cef5090cSJed Brown {
3186cef5090cSJed Brown   PetscFunctionBegin;
3187cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3188cef5090cSJed Brown   PetscValidIntPointer(fails,2);
3189cef5090cSJed Brown   *fails = ts->num_snes_failures;
3190cef5090cSJed Brown   PetscFunctionReturn(0);
3191cef5090cSJed Brown }
3192cef5090cSJed Brown 
3193cef5090cSJed Brown #undef __FUNCT__
3194cef5090cSJed Brown #define __FUNCT__ "TSSetMaxStepRejections"
3195cef5090cSJed Brown /*@
3196cef5090cSJed Brown    TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails
3197cef5090cSJed Brown 
3198cef5090cSJed Brown    Not Collective
3199cef5090cSJed Brown 
3200cef5090cSJed Brown    Input Parameter:
3201cef5090cSJed Brown +  ts - TS context
3202cef5090cSJed Brown -  rejects - maximum number of rejected steps, pass -1 for unlimited
3203cef5090cSJed Brown 
3204cef5090cSJed Brown    Notes:
3205cef5090cSJed Brown    The counter is reset to zero for each step
3206cef5090cSJed Brown 
3207cef5090cSJed Brown    Options Database Key:
3208cef5090cSJed Brown  .  -ts_max_reject - Maximum number of step rejections before a step fails
3209cef5090cSJed Brown 
3210cef5090cSJed Brown    Level: intermediate
3211cef5090cSJed Brown 
3212cef5090cSJed Brown .keywords: TS, set, maximum, number
3213cef5090cSJed Brown 
32145ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3215cef5090cSJed Brown @*/
3216cef5090cSJed Brown PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
3217cef5090cSJed Brown {
3218cef5090cSJed Brown   PetscFunctionBegin;
3219cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3220cef5090cSJed Brown   ts->max_reject = rejects;
3221cef5090cSJed Brown   PetscFunctionReturn(0);
3222cef5090cSJed Brown }
3223cef5090cSJed Brown 
3224cef5090cSJed Brown #undef __FUNCT__
3225cef5090cSJed Brown #define __FUNCT__ "TSSetMaxSNESFailures"
3226cef5090cSJed Brown /*@
3227cef5090cSJed Brown    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves
3228cef5090cSJed Brown 
3229cef5090cSJed Brown    Not Collective
3230cef5090cSJed Brown 
3231cef5090cSJed Brown    Input Parameter:
3232cef5090cSJed Brown +  ts - TS context
3233cef5090cSJed Brown -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited
3234cef5090cSJed Brown 
3235cef5090cSJed Brown    Notes:
3236cef5090cSJed Brown    The counter is reset to zero for each successive call to TSSolve().
3237cef5090cSJed Brown 
3238cef5090cSJed Brown    Options Database Key:
3239cef5090cSJed Brown  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures
3240cef5090cSJed Brown 
3241cef5090cSJed Brown    Level: intermediate
3242cef5090cSJed Brown 
3243cef5090cSJed Brown .keywords: TS, set, maximum, number
3244cef5090cSJed Brown 
32455ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
3246cef5090cSJed Brown @*/
3247cef5090cSJed Brown PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
3248cef5090cSJed Brown {
3249cef5090cSJed Brown   PetscFunctionBegin;
3250cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3251cef5090cSJed Brown   ts->max_snes_failures = fails;
3252cef5090cSJed Brown   PetscFunctionReturn(0);
3253cef5090cSJed Brown }
3254cef5090cSJed Brown 
3255cef5090cSJed Brown #undef __FUNCT__
3256cef5090cSJed Brown #define __FUNCT__ "TSSetErrorIfStepFails()"
3257cef5090cSJed Brown /*@
3258cef5090cSJed Brown    TSSetErrorIfStepFails - Error if no step succeeds
3259cef5090cSJed Brown 
3260cef5090cSJed Brown    Not Collective
3261cef5090cSJed Brown 
3262cef5090cSJed Brown    Input Parameter:
3263cef5090cSJed Brown +  ts - TS context
3264cef5090cSJed Brown -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure
3265cef5090cSJed Brown 
3266cef5090cSJed Brown    Options Database Key:
3267cef5090cSJed Brown  .  -ts_error_if_step_fails - Error if no step succeeds
3268cef5090cSJed Brown 
3269cef5090cSJed Brown    Level: intermediate
3270cef5090cSJed Brown 
3271cef5090cSJed Brown .keywords: TS, set, error
3272cef5090cSJed Brown 
32735ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3274cef5090cSJed Brown @*/
3275cef5090cSJed Brown PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
3276cef5090cSJed Brown {
3277cef5090cSJed Brown   PetscFunctionBegin;
3278cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3279cef5090cSJed Brown   ts->errorifstepfailed = err;
3280cef5090cSJed Brown   PetscFunctionReturn(0);
3281cef5090cSJed Brown }
3282cef5090cSJed Brown 
3283cef5090cSJed Brown #undef __FUNCT__
3284fb1732b5SBarry Smith #define __FUNCT__ "TSMonitorSolutionBinary"
3285fb1732b5SBarry Smith /*@C
3286fb1732b5SBarry Smith    TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file
3287fb1732b5SBarry Smith 
3288fb1732b5SBarry Smith    Collective on TS
3289fb1732b5SBarry Smith 
3290fb1732b5SBarry Smith    Input Parameters:
3291fb1732b5SBarry Smith +  ts - the TS context
3292fb1732b5SBarry Smith .  step - current time-step
3293fb1732b5SBarry Smith .  ptime - current time
3294ed81e22dSJed Brown .  x - current state
3295fb1732b5SBarry Smith -  viewer - binary viewer
3296fb1732b5SBarry Smith 
3297fb1732b5SBarry Smith    Level: intermediate
3298fb1732b5SBarry Smith 
3299fb1732b5SBarry Smith .keywords: TS,  vector, monitor, view
3300fb1732b5SBarry Smith 
3301fb1732b5SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3302fb1732b5SBarry Smith @*/
3303ed81e22dSJed Brown PetscErrorCode  TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer)
3304fb1732b5SBarry Smith {
3305fb1732b5SBarry Smith   PetscErrorCode       ierr;
3306ed81e22dSJed Brown   PetscViewer          v = (PetscViewer)viewer;
3307fb1732b5SBarry Smith 
3308fb1732b5SBarry Smith   PetscFunctionBegin;
3309ed81e22dSJed Brown   ierr = VecView(x,v);CHKERRQ(ierr);
3310ed81e22dSJed Brown   PetscFunctionReturn(0);
3311ed81e22dSJed Brown }
3312ed81e22dSJed Brown 
3313ed81e22dSJed Brown #undef __FUNCT__
3314ed81e22dSJed Brown #define __FUNCT__ "TSMonitorSolutionVTK"
3315ed81e22dSJed Brown /*@C
3316ed81e22dSJed Brown    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.
3317ed81e22dSJed Brown 
3318ed81e22dSJed Brown    Collective on TS
3319ed81e22dSJed Brown 
3320ed81e22dSJed Brown    Input Parameters:
3321ed81e22dSJed Brown +  ts - the TS context
3322ed81e22dSJed Brown .  step - current time-step
3323ed81e22dSJed Brown .  ptime - current time
3324ed81e22dSJed Brown .  x - current state
3325ed81e22dSJed Brown -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3326ed81e22dSJed Brown 
3327ed81e22dSJed Brown    Level: intermediate
3328ed81e22dSJed Brown 
3329ed81e22dSJed Brown    Notes:
3330ed81e22dSJed 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.
3331ed81e22dSJed Brown    These are named according to the file name template.
3332ed81e22dSJed Brown 
3333ed81e22dSJed Brown    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().
3334ed81e22dSJed Brown 
3335ed81e22dSJed Brown .keywords: TS,  vector, monitor, view
3336ed81e22dSJed Brown 
3337ed81e22dSJed Brown .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3338ed81e22dSJed Brown @*/
3339ed81e22dSJed Brown PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate)
3340ed81e22dSJed Brown {
3341ed81e22dSJed Brown   PetscErrorCode ierr;
3342ed81e22dSJed Brown   char           filename[PETSC_MAX_PATH_LEN];
3343ed81e22dSJed Brown   PetscViewer    viewer;
3344ed81e22dSJed Brown 
3345ed81e22dSJed Brown   PetscFunctionBegin;
33468caf3d72SBarry Smith   ierr = PetscSNPrintf(filename,sizeof(filename),(const char*)filenametemplate,step);CHKERRQ(ierr);
3347ed81e22dSJed Brown   ierr = PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
3348fb1732b5SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
3349ed81e22dSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3350ed81e22dSJed Brown   PetscFunctionReturn(0);
3351ed81e22dSJed Brown }
3352ed81e22dSJed Brown 
3353ed81e22dSJed Brown #undef __FUNCT__
3354ed81e22dSJed Brown #define __FUNCT__ "TSMonitorSolutionVTKDestroy"
3355ed81e22dSJed Brown /*@C
3356ed81e22dSJed Brown    TSMonitorSolutionVTKDestroy - Destroy context for monitoring
3357ed81e22dSJed Brown 
3358ed81e22dSJed Brown    Collective on TS
3359ed81e22dSJed Brown 
3360ed81e22dSJed Brown    Input Parameters:
3361ed81e22dSJed Brown .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3362ed81e22dSJed Brown 
3363ed81e22dSJed Brown    Level: intermediate
3364ed81e22dSJed Brown 
3365ed81e22dSJed Brown    Note:
3366ed81e22dSJed Brown    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().
3367ed81e22dSJed Brown 
3368ed81e22dSJed Brown .keywords: TS,  vector, monitor, view
3369ed81e22dSJed Brown 
3370ed81e22dSJed Brown .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
3371ed81e22dSJed Brown @*/
3372ed81e22dSJed Brown PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
3373ed81e22dSJed Brown {
3374ed81e22dSJed Brown   PetscErrorCode ierr;
3375ed81e22dSJed Brown 
3376ed81e22dSJed Brown   PetscFunctionBegin;
3377ed81e22dSJed Brown   ierr = PetscFree(*(char**)filenametemplate);CHKERRQ(ierr);
3378fb1732b5SBarry Smith   PetscFunctionReturn(0);
3379fb1732b5SBarry Smith }
3380fb1732b5SBarry Smith 
338184df9cb4SJed Brown #undef __FUNCT__
338284df9cb4SJed Brown #define __FUNCT__ "TSGetAdapt"
338384df9cb4SJed Brown /*@
338484df9cb4SJed Brown    TSGetAdapt - Get the adaptive controller context for the current method
338584df9cb4SJed Brown 
3386ed81e22dSJed Brown    Collective on TS if controller has not been created yet
338784df9cb4SJed Brown 
338884df9cb4SJed Brown    Input Arguments:
3389ed81e22dSJed Brown .  ts - time stepping context
339084df9cb4SJed Brown 
339184df9cb4SJed Brown    Output Arguments:
3392ed81e22dSJed Brown .  adapt - adaptive controller
339384df9cb4SJed Brown 
339484df9cb4SJed Brown    Level: intermediate
339584df9cb4SJed Brown 
3396ed81e22dSJed Brown .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
339784df9cb4SJed Brown @*/
339884df9cb4SJed Brown PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
339984df9cb4SJed Brown {
340084df9cb4SJed Brown   PetscErrorCode ierr;
340184df9cb4SJed Brown 
340284df9cb4SJed Brown   PetscFunctionBegin;
340384df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
340484df9cb4SJed Brown   PetscValidPointer(adapt,2);
340584df9cb4SJed Brown   if (!ts->adapt) {
340684df9cb4SJed Brown     ierr = TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);CHKERRQ(ierr);
34071c3436cfSJed Brown     ierr = PetscLogObjectParent(ts,ts->adapt);CHKERRQ(ierr);
34081c3436cfSJed Brown     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);CHKERRQ(ierr);
340984df9cb4SJed Brown   }
341084df9cb4SJed Brown   *adapt = ts->adapt;
341184df9cb4SJed Brown   PetscFunctionReturn(0);
341284df9cb4SJed Brown }
3413d6ebe24aSShri Abhyankar 
3414d6ebe24aSShri Abhyankar #undef __FUNCT__
34151c3436cfSJed Brown #define __FUNCT__ "TSSetTolerances"
34161c3436cfSJed Brown /*@
34171c3436cfSJed Brown    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller
34181c3436cfSJed Brown 
34191c3436cfSJed Brown    Logically Collective
34201c3436cfSJed Brown 
34211c3436cfSJed Brown    Input Arguments:
34221c3436cfSJed Brown +  ts - time integration context
34231c3436cfSJed Brown .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
34241c3436cfSJed Brown .  vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present
34251c3436cfSJed Brown .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
34261c3436cfSJed Brown -  vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present
34271c3436cfSJed Brown 
34281c3436cfSJed Brown    Level: beginner
34291c3436cfSJed Brown 
3430c5033834SJed Brown .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
34311c3436cfSJed Brown @*/
34321c3436cfSJed Brown PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
34331c3436cfSJed Brown {
34341c3436cfSJed Brown   PetscErrorCode ierr;
34351c3436cfSJed Brown 
34361c3436cfSJed Brown   PetscFunctionBegin;
3437c5033834SJed Brown   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
34381c3436cfSJed Brown   if (vatol) {
34391c3436cfSJed Brown     ierr = PetscObjectReference((PetscObject)vatol);CHKERRQ(ierr);
34401c3436cfSJed Brown     ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
34411c3436cfSJed Brown     ts->vatol = vatol;
34421c3436cfSJed Brown   }
3443c5033834SJed Brown   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
34441c3436cfSJed Brown   if (vrtol) {
34451c3436cfSJed Brown     ierr = PetscObjectReference((PetscObject)vrtol);CHKERRQ(ierr);
34461c3436cfSJed Brown     ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
34471c3436cfSJed Brown     ts->vrtol = vrtol;
34481c3436cfSJed Brown   }
34491c3436cfSJed Brown   PetscFunctionReturn(0);
34501c3436cfSJed Brown }
34511c3436cfSJed Brown 
34521c3436cfSJed Brown #undef __FUNCT__
3453c5033834SJed Brown #define __FUNCT__ "TSGetTolerances"
3454c5033834SJed Brown /*@
3455c5033834SJed Brown    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
3456c5033834SJed Brown 
3457c5033834SJed Brown    Logically Collective
3458c5033834SJed Brown 
3459c5033834SJed Brown    Input Arguments:
3460c5033834SJed Brown .  ts - time integration context
3461c5033834SJed Brown 
3462c5033834SJed Brown    Output Arguments:
3463c5033834SJed Brown +  atol - scalar absolute tolerances, PETSC_NULL to ignore
3464c5033834SJed Brown .  vatol - vector of absolute tolerances, PETSC_NULL to ignore
3465c5033834SJed Brown .  rtol - scalar relative tolerances, PETSC_NULL to ignore
3466c5033834SJed Brown -  vrtol - vector of relative tolerances, PETSC_NULL to ignore
3467c5033834SJed Brown 
3468c5033834SJed Brown    Level: beginner
3469c5033834SJed Brown 
3470c5033834SJed Brown .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
3471c5033834SJed Brown @*/
3472c5033834SJed Brown PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
3473c5033834SJed Brown {
3474c5033834SJed Brown 
3475c5033834SJed Brown   PetscFunctionBegin;
3476c5033834SJed Brown   if (atol)  *atol  = ts->atol;
3477c5033834SJed Brown   if (vatol) *vatol = ts->vatol;
3478c5033834SJed Brown   if (rtol)  *rtol  = ts->rtol;
3479c5033834SJed Brown   if (vrtol) *vrtol = ts->vrtol;
3480c5033834SJed Brown   PetscFunctionReturn(0);
3481c5033834SJed Brown }
3482c5033834SJed Brown 
3483c5033834SJed Brown #undef __FUNCT__
34841c3436cfSJed Brown #define __FUNCT__ "TSErrorNormWRMS"
34851c3436cfSJed Brown /*@
34861c3436cfSJed Brown    TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state
34871c3436cfSJed Brown 
34881c3436cfSJed Brown    Collective on TS
34891c3436cfSJed Brown 
34901c3436cfSJed Brown    Input Arguments:
34911c3436cfSJed Brown +  ts - time stepping context
34921c3436cfSJed Brown -  Y - state vector to be compared to ts->vec_sol
34931c3436cfSJed Brown 
34941c3436cfSJed Brown    Output Arguments:
34951c3436cfSJed Brown .  norm - weighted norm, a value of 1.0 is considered small
34961c3436cfSJed Brown 
34971c3436cfSJed Brown    Level: developer
34981c3436cfSJed Brown 
34991c3436cfSJed Brown .seealso: TSSetTolerances()
35001c3436cfSJed Brown @*/
35011c3436cfSJed Brown PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
35021c3436cfSJed Brown {
35031c3436cfSJed Brown   PetscErrorCode ierr;
35041c3436cfSJed Brown   PetscInt i,n,N;
35051c3436cfSJed Brown   const PetscScalar *x,*y;
35061c3436cfSJed Brown   Vec X;
35071c3436cfSJed Brown   PetscReal sum,gsum;
35081c3436cfSJed Brown 
35091c3436cfSJed Brown   PetscFunctionBegin;
35101c3436cfSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
35111c3436cfSJed Brown   PetscValidHeaderSpecific(Y,VEC_CLASSID,2);
35121c3436cfSJed Brown   PetscValidPointer(norm,3);
35131c3436cfSJed Brown   X = ts->vec_sol;
35141c3436cfSJed Brown   PetscCheckSameTypeAndComm(X,1,Y,2);
35151c3436cfSJed Brown   if (X == Y) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector");
35161c3436cfSJed Brown 
35171c3436cfSJed Brown   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
35181c3436cfSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
35191c3436cfSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
35201c3436cfSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
35211c3436cfSJed Brown   sum = 0.;
3522ceefc113SJed Brown   if (ts->vatol && ts->vrtol) {
3523ceefc113SJed Brown     const PetscScalar *atol,*rtol;
3524ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3525ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3526ceefc113SJed Brown     for (i=0; i<n; i++) {
3527ceefc113SJed Brown       PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3528ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3529ceefc113SJed Brown     }
3530ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3531ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3532ceefc113SJed Brown   } else if (ts->vatol) {       /* vector atol, scalar rtol */
3533ceefc113SJed Brown     const PetscScalar *atol;
3534ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3535ceefc113SJed Brown     for (i=0; i<n; i++) {
3536ceefc113SJed Brown       PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3537ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3538ceefc113SJed Brown     }
3539ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3540ceefc113SJed Brown   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
3541ceefc113SJed Brown     const PetscScalar *rtol;
3542ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3543ceefc113SJed Brown     for (i=0; i<n; i++) {
3544ceefc113SJed Brown       PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3545ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3546ceefc113SJed Brown     }
3547ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3548ceefc113SJed Brown   } else {                      /* scalar atol, scalar rtol */
35491c3436cfSJed Brown     for (i=0; i<n; i++) {
35501c3436cfSJed Brown       PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
35511c3436cfSJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
35521c3436cfSJed Brown     }
3553ceefc113SJed Brown   }
35541c3436cfSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
35551c3436cfSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
35561c3436cfSJed Brown 
35571c3436cfSJed Brown   ierr = MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);CHKERRQ(ierr);
35581c3436cfSJed Brown   *norm = PetscSqrtReal(gsum / N);
35591c3436cfSJed Brown   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
35601c3436cfSJed Brown   PetscFunctionReturn(0);
35611c3436cfSJed Brown }
35621c3436cfSJed Brown 
35631c3436cfSJed Brown #undef __FUNCT__
35648d59e960SJed Brown #define __FUNCT__ "TSSetCFLTimeLocal"
35658d59e960SJed Brown /*@
35668d59e960SJed Brown    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
35678d59e960SJed Brown 
35688d59e960SJed Brown    Logically Collective on TS
35698d59e960SJed Brown 
35708d59e960SJed Brown    Input Arguments:
35718d59e960SJed Brown +  ts - time stepping context
35728d59e960SJed Brown -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)
35738d59e960SJed Brown 
35748d59e960SJed Brown    Note:
35758d59e960SJed Brown    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
35768d59e960SJed Brown 
35778d59e960SJed Brown    Level: intermediate
35788d59e960SJed Brown 
35798d59e960SJed Brown .seealso: TSGetCFLTime(), TSADAPTCFL
35808d59e960SJed Brown @*/
35818d59e960SJed Brown PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
35828d59e960SJed Brown {
35838d59e960SJed Brown 
35848d59e960SJed Brown   PetscFunctionBegin;
35858d59e960SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
35868d59e960SJed Brown   ts->cfltime_local = cfltime;
35878d59e960SJed Brown   ts->cfltime = -1.;
35888d59e960SJed Brown   PetscFunctionReturn(0);
35898d59e960SJed Brown }
35908d59e960SJed Brown 
35918d59e960SJed Brown #undef __FUNCT__
35928d59e960SJed Brown #define __FUNCT__ "TSGetCFLTime"
35938d59e960SJed Brown /*@
35948d59e960SJed Brown    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
35958d59e960SJed Brown 
35968d59e960SJed Brown    Collective on TS
35978d59e960SJed Brown 
35988d59e960SJed Brown    Input Arguments:
35998d59e960SJed Brown .  ts - time stepping context
36008d59e960SJed Brown 
36018d59e960SJed Brown    Output Arguments:
36028d59e960SJed Brown .  cfltime - maximum stable time step for forward Euler
36038d59e960SJed Brown 
36048d59e960SJed Brown    Level: advanced
36058d59e960SJed Brown 
36068d59e960SJed Brown .seealso: TSSetCFLTimeLocal()
36078d59e960SJed Brown @*/
36088d59e960SJed Brown PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
36098d59e960SJed Brown {
36108d59e960SJed Brown   PetscErrorCode ierr;
36118d59e960SJed Brown 
36128d59e960SJed Brown   PetscFunctionBegin;
36138d59e960SJed Brown   if (ts->cfltime < 0) {
36148d59e960SJed Brown     ierr = MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr);
36158d59e960SJed Brown   }
36168d59e960SJed Brown   *cfltime = ts->cfltime;
36178d59e960SJed Brown   PetscFunctionReturn(0);
36188d59e960SJed Brown }
36198d59e960SJed Brown 
36208d59e960SJed Brown #undef __FUNCT__
3621d6ebe24aSShri Abhyankar #define __FUNCT__ "TSVISetVariableBounds"
3622d6ebe24aSShri Abhyankar /*@
3623d6ebe24aSShri Abhyankar    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
3624d6ebe24aSShri Abhyankar 
3625d6ebe24aSShri Abhyankar    Input Parameters:
3626d6ebe24aSShri Abhyankar .  ts   - the TS context.
3627d6ebe24aSShri Abhyankar .  xl   - lower bound.
3628d6ebe24aSShri Abhyankar .  xu   - upper bound.
3629d6ebe24aSShri Abhyankar 
3630d6ebe24aSShri Abhyankar    Notes:
3631d6ebe24aSShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
363233c0c2ddSJed Brown    SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp().
3633d6ebe24aSShri Abhyankar 
36342bd2b0e6SSatish Balay    Level: advanced
36352bd2b0e6SSatish Balay 
3636d6ebe24aSShri Abhyankar @*/
3637d6ebe24aSShri Abhyankar PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
3638d6ebe24aSShri Abhyankar {
3639d6ebe24aSShri Abhyankar   PetscErrorCode ierr;
3640d6ebe24aSShri Abhyankar   SNES           snes;
3641d6ebe24aSShri Abhyankar 
3642d6ebe24aSShri Abhyankar   PetscFunctionBegin;
3643d6ebe24aSShri Abhyankar   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3644d6ebe24aSShri Abhyankar   ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr);
3645d6ebe24aSShri Abhyankar   PetscFunctionReturn(0);
3646d6ebe24aSShri Abhyankar }
3647d6ebe24aSShri Abhyankar 
3648325fc9f4SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
3649c6db04a5SJed Brown #include <mex.h>
3650325fc9f4SBarry Smith 
3651325fc9f4SBarry Smith typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
3652325fc9f4SBarry Smith 
3653325fc9f4SBarry Smith #undef __FUNCT__
3654325fc9f4SBarry Smith #define __FUNCT__ "TSComputeFunction_Matlab"
3655325fc9f4SBarry Smith /*
3656325fc9f4SBarry Smith    TSComputeFunction_Matlab - Calls the function that has been set with
3657325fc9f4SBarry Smith                          TSSetFunctionMatlab().
3658325fc9f4SBarry Smith 
3659325fc9f4SBarry Smith    Collective on TS
3660325fc9f4SBarry Smith 
3661325fc9f4SBarry Smith    Input Parameters:
3662325fc9f4SBarry Smith +  snes - the TS context
3663325fc9f4SBarry Smith -  x - input vector
3664325fc9f4SBarry Smith 
3665325fc9f4SBarry Smith    Output Parameter:
3666325fc9f4SBarry Smith .  y - function vector, as set by TSSetFunction()
3667325fc9f4SBarry Smith 
3668325fc9f4SBarry Smith    Notes:
3669325fc9f4SBarry Smith    TSComputeFunction() is typically used within nonlinear solvers
3670325fc9f4SBarry Smith    implementations, so most users would not generally call this routine
3671325fc9f4SBarry Smith    themselves.
3672325fc9f4SBarry Smith 
3673325fc9f4SBarry Smith    Level: developer
3674325fc9f4SBarry Smith 
3675325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
3676325fc9f4SBarry Smith 
3677325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3678325fc9f4SBarry Smith */
36797087cfbeSBarry Smith PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
3680325fc9f4SBarry Smith {
3681325fc9f4SBarry Smith   PetscErrorCode   ierr;
3682325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3683325fc9f4SBarry Smith   int              nlhs = 1,nrhs = 7;
3684325fc9f4SBarry Smith   mxArray          *plhs[1],*prhs[7];
3685325fc9f4SBarry Smith   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
3686325fc9f4SBarry Smith 
3687325fc9f4SBarry Smith   PetscFunctionBegin;
3688325fc9f4SBarry Smith   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
3689325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3690325fc9f4SBarry Smith   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
3691325fc9f4SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
3692325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,x,3);
3693325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,y,5);
3694325fc9f4SBarry Smith 
3695325fc9f4SBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3696325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3697880f3077SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
3698325fc9f4SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
3699325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3700325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar(time);
3701325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
3702325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3703325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)ly);
3704325fc9f4SBarry Smith   prhs[5] =  mxCreateString(sctx->funcname);
3705325fc9f4SBarry Smith   prhs[6] =  sctx->ctx;
3706325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
3707325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3708325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
3709325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
3710325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
3711325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
3712325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
3713325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
3714325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
3715325fc9f4SBarry Smith   PetscFunctionReturn(0);
3716325fc9f4SBarry Smith }
3717325fc9f4SBarry Smith 
3718325fc9f4SBarry Smith 
3719325fc9f4SBarry Smith #undef __FUNCT__
3720325fc9f4SBarry Smith #define __FUNCT__ "TSSetFunctionMatlab"
3721325fc9f4SBarry Smith /*
3722325fc9f4SBarry Smith    TSSetFunctionMatlab - Sets the function evaluation routine and function
3723325fc9f4SBarry Smith    vector for use by the TS routines in solving ODEs
3724e3c5b3baSBarry Smith    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
3725325fc9f4SBarry Smith 
3726325fc9f4SBarry Smith    Logically Collective on TS
3727325fc9f4SBarry Smith 
3728325fc9f4SBarry Smith    Input Parameters:
3729325fc9f4SBarry Smith +  ts - the TS context
3730325fc9f4SBarry Smith -  func - function evaluation routine
3731325fc9f4SBarry Smith 
3732325fc9f4SBarry Smith    Calling sequence of func:
3733325fc9f4SBarry Smith $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
3734325fc9f4SBarry Smith 
3735325fc9f4SBarry Smith    Level: beginner
3736325fc9f4SBarry Smith 
3737325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
3738325fc9f4SBarry Smith 
3739325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3740325fc9f4SBarry Smith */
3741cdcf91faSSean Farley PetscErrorCode  TSSetFunctionMatlab(TS ts,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 = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
3757325fc9f4SBarry Smith   PetscFunctionReturn(0);
3758325fc9f4SBarry Smith }
3759325fc9f4SBarry Smith 
3760325fc9f4SBarry Smith #undef __FUNCT__
3761325fc9f4SBarry Smith #define __FUNCT__ "TSComputeJacobian_Matlab"
3762325fc9f4SBarry Smith /*
3763325fc9f4SBarry Smith    TSComputeJacobian_Matlab - Calls the function that has been set with
3764325fc9f4SBarry Smith                          TSSetJacobianMatlab().
3765325fc9f4SBarry Smith 
3766325fc9f4SBarry Smith    Collective on TS
3767325fc9f4SBarry Smith 
3768325fc9f4SBarry Smith    Input Parameters:
3769cdcf91faSSean Farley +  ts - the TS context
3770325fc9f4SBarry Smith .  x - input vector
3771325fc9f4SBarry Smith .  A, B - the matrices
3772325fc9f4SBarry Smith -  ctx - user context
3773325fc9f4SBarry Smith 
3774325fc9f4SBarry Smith    Output Parameter:
3775325fc9f4SBarry Smith .  flag - structure of the matrix
3776325fc9f4SBarry Smith 
3777325fc9f4SBarry Smith    Level: developer
3778325fc9f4SBarry Smith 
3779325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
3780325fc9f4SBarry Smith 
3781325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3782325fc9f4SBarry Smith @*/
3783cdcf91faSSean Farley PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
3784325fc9f4SBarry Smith {
3785325fc9f4SBarry Smith   PetscErrorCode  ierr;
3786325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3787325fc9f4SBarry Smith   int             nlhs = 2,nrhs = 9;
3788325fc9f4SBarry Smith   mxArray         *plhs[2],*prhs[9];
3789325fc9f4SBarry Smith   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
3790325fc9f4SBarry Smith 
3791325fc9f4SBarry Smith   PetscFunctionBegin;
3792cdcf91faSSean Farley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3793325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3794325fc9f4SBarry Smith 
3795325fc9f4SBarry Smith   /* call Matlab function in ctx with arguments x and y */
3796325fc9f4SBarry Smith 
3797cdcf91faSSean Farley   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
3798325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3799325fc9f4SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
3800325fc9f4SBarry Smith   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
3801325fc9f4SBarry Smith   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
3802325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3803325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)time);
3804325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
3805325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3806325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)shift);
3807325fc9f4SBarry Smith   prhs[5] =  mxCreateDoubleScalar((double)lA);
3808325fc9f4SBarry Smith   prhs[6] =  mxCreateDoubleScalar((double)lB);
3809325fc9f4SBarry Smith   prhs[7] =  mxCreateString(sctx->funcname);
3810325fc9f4SBarry Smith   prhs[8] =  sctx->ctx;
3811325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
3812325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3813325fc9f4SBarry Smith   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
3814325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
3815325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
3816325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
3817325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
3818325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
3819325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
3820325fc9f4SBarry Smith   mxDestroyArray(prhs[6]);
3821325fc9f4SBarry Smith   mxDestroyArray(prhs[7]);
3822325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
3823325fc9f4SBarry Smith   mxDestroyArray(plhs[1]);
3824325fc9f4SBarry Smith   PetscFunctionReturn(0);
3825325fc9f4SBarry Smith }
3826325fc9f4SBarry Smith 
3827325fc9f4SBarry Smith 
3828325fc9f4SBarry Smith #undef __FUNCT__
3829325fc9f4SBarry Smith #define __FUNCT__ "TSSetJacobianMatlab"
3830325fc9f4SBarry Smith /*
3831325fc9f4SBarry Smith    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
3832e3c5b3baSBarry 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
3833325fc9f4SBarry Smith 
3834325fc9f4SBarry Smith    Logically Collective on TS
3835325fc9f4SBarry Smith 
3836325fc9f4SBarry Smith    Input Parameters:
3837cdcf91faSSean Farley +  ts - the TS context
3838325fc9f4SBarry Smith .  A,B - Jacobian matrices
3839325fc9f4SBarry Smith .  func - function evaluation routine
3840325fc9f4SBarry Smith -  ctx - user context
3841325fc9f4SBarry Smith 
3842325fc9f4SBarry Smith    Calling sequence of func:
3843cdcf91faSSean Farley $    flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
3844325fc9f4SBarry Smith 
3845325fc9f4SBarry Smith 
3846325fc9f4SBarry Smith    Level: developer
3847325fc9f4SBarry Smith 
3848325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
3849325fc9f4SBarry Smith 
3850325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3851325fc9f4SBarry Smith */
3852cdcf91faSSean Farley PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
3853325fc9f4SBarry Smith {
3854325fc9f4SBarry Smith   PetscErrorCode    ierr;
3855325fc9f4SBarry Smith   TSMatlabContext *sctx;
3856325fc9f4SBarry Smith 
3857325fc9f4SBarry Smith   PetscFunctionBegin;
3858325fc9f4SBarry Smith   /* currently sctx is memory bleed */
3859325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3860325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3861325fc9f4SBarry Smith   /*
3862325fc9f4SBarry Smith      This should work, but it doesn't
3863325fc9f4SBarry Smith   sctx->ctx = ctx;
3864325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
3865325fc9f4SBarry Smith   */
3866325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
3867cdcf91faSSean Farley   ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
3868325fc9f4SBarry Smith   PetscFunctionReturn(0);
3869325fc9f4SBarry Smith }
3870325fc9f4SBarry Smith 
3871b5b1a830SBarry Smith #undef __FUNCT__
3872b5b1a830SBarry Smith #define __FUNCT__ "TSMonitor_Matlab"
3873b5b1a830SBarry Smith /*
3874b5b1a830SBarry Smith    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
3875b5b1a830SBarry Smith 
3876b5b1a830SBarry Smith    Collective on TS
3877b5b1a830SBarry Smith 
3878b5b1a830SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3879b5b1a830SBarry Smith @*/
3880cdcf91faSSean Farley PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
3881b5b1a830SBarry Smith {
3882b5b1a830SBarry Smith   PetscErrorCode  ierr;
3883b5b1a830SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3884a530c242SBarry Smith   int             nlhs = 1,nrhs = 6;
3885b5b1a830SBarry Smith   mxArray         *plhs[1],*prhs[6];
3886b5b1a830SBarry Smith   long long int   lx = 0,ls = 0;
3887b5b1a830SBarry Smith 
3888b5b1a830SBarry Smith   PetscFunctionBegin;
3889cdcf91faSSean Farley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3890b5b1a830SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
3891b5b1a830SBarry Smith 
3892cdcf91faSSean Farley   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
3893b5b1a830SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3894b5b1a830SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3895b5b1a830SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)it);
3896b5b1a830SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)time);
3897b5b1a830SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lx);
3898b5b1a830SBarry Smith   prhs[4] =  mxCreateString(sctx->funcname);
3899b5b1a830SBarry Smith   prhs[5] =  sctx->ctx;
3900b5b1a830SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
3901b5b1a830SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3902b5b1a830SBarry Smith   mxDestroyArray(prhs[0]);
3903b5b1a830SBarry Smith   mxDestroyArray(prhs[1]);
3904b5b1a830SBarry Smith   mxDestroyArray(prhs[2]);
3905b5b1a830SBarry Smith   mxDestroyArray(prhs[3]);
3906b5b1a830SBarry Smith   mxDestroyArray(prhs[4]);
3907b5b1a830SBarry Smith   mxDestroyArray(plhs[0]);
3908b5b1a830SBarry Smith   PetscFunctionReturn(0);
3909b5b1a830SBarry Smith }
3910b5b1a830SBarry Smith 
3911b5b1a830SBarry Smith 
3912b5b1a830SBarry Smith #undef __FUNCT__
3913b5b1a830SBarry Smith #define __FUNCT__ "TSMonitorSetMatlab"
3914b5b1a830SBarry Smith /*
3915b5b1a830SBarry Smith    TSMonitorSetMatlab - Sets the monitor function from Matlab
3916b5b1a830SBarry Smith 
3917b5b1a830SBarry Smith    Level: developer
3918b5b1a830SBarry Smith 
3919b5b1a830SBarry Smith .keywords: TS, nonlinear, set, function
3920b5b1a830SBarry Smith 
3921b5b1a830SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3922b5b1a830SBarry Smith */
3923cdcf91faSSean Farley PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
3924b5b1a830SBarry Smith {
3925b5b1a830SBarry Smith   PetscErrorCode    ierr;
3926b5b1a830SBarry Smith   TSMatlabContext *sctx;
3927b5b1a830SBarry Smith 
3928b5b1a830SBarry Smith   PetscFunctionBegin;
3929b5b1a830SBarry Smith   /* currently sctx is memory bleed */
3930b5b1a830SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3931b5b1a830SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3932b5b1a830SBarry Smith   /*
3933b5b1a830SBarry Smith      This should work, but it doesn't
3934b5b1a830SBarry Smith   sctx->ctx = ctx;
3935b5b1a830SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
3936b5b1a830SBarry Smith   */
3937b5b1a830SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
3938cdcf91faSSean Farley   ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
3939b5b1a830SBarry Smith   PetscFunctionReturn(0);
3940b5b1a830SBarry Smith }
3941325fc9f4SBarry Smith #endif
3942b3603a34SBarry Smith 
3943b3603a34SBarry Smith #undef __FUNCT__
3944b3603a34SBarry Smith #define __FUNCT__ "TSMonitorSolutionODE"
3945b3603a34SBarry Smith /*@C
3946b3603a34SBarry Smith    TSMonitorSolutionODE - Monitors progress of the TS solvers by plotting each component of the solution vector
3947b3603a34SBarry Smith        in a time based line graph
3948b3603a34SBarry Smith 
3949b3603a34SBarry Smith    Collective on TS
3950b3603a34SBarry Smith 
3951b3603a34SBarry Smith    Input Parameters:
3952b3603a34SBarry Smith +  ts - the TS context
3953b3603a34SBarry Smith .  step - current time-step
3954b3603a34SBarry Smith .  ptime - current time
3955b3603a34SBarry Smith -  lg - a line graph object
3956b3603a34SBarry Smith 
3957b3603a34SBarry Smith    Level: intermediate
3958b3603a34SBarry Smith 
3959b3603a34SBarry Smith     Notes: only for sequential solves
3960b3603a34SBarry Smith 
3961b3603a34SBarry Smith .keywords: TS,  vector, monitor, view
3962b3603a34SBarry Smith 
3963b3603a34SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3964b3603a34SBarry Smith @*/
3965b3603a34SBarry Smith PetscErrorCode  TSMonitorSolutionODE(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
3966b3603a34SBarry Smith {
3967b3603a34SBarry Smith   PetscErrorCode    ierr;
3968b3603a34SBarry Smith   PetscDrawLG       lg = (PetscDrawLG)dummy;
3969b3603a34SBarry Smith   const PetscScalar *yy;
3970b3603a34SBarry Smith 
3971b3603a34SBarry Smith   PetscFunctionBegin;
3972b3603a34SBarry Smith   if (!step) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
3973b3603a34SBarry Smith   ierr = VecGetArrayRead(x,&yy);CHKERRQ(ierr);
3974aa39b21eSBarry Smith   ierr = PetscDrawLGAddCommonPoint(lg,ptime,yy);CHKERRQ(ierr);
3975b3603a34SBarry Smith   ierr = VecRestoreArrayRead(x,&yy);CHKERRQ(ierr);
3976b3603a34SBarry Smith   ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
3977b3603a34SBarry Smith   PetscFunctionReturn(0);
3978b3603a34SBarry Smith }
3979b3603a34SBarry Smith 
3980b3603a34SBarry Smith #undef __FUNCT__
3981b3603a34SBarry Smith #define __FUNCT__ "TSMonitorSolutionODEDestroy"
3982b3603a34SBarry Smith /*@C
3983b3603a34SBarry Smith    TSMonitorSolutionODEDestroy - Destroys the monitor context for TSMonitorSolutionODE()
3984b3603a34SBarry Smith 
3985b3603a34SBarry Smith    Collective on TS
3986b3603a34SBarry Smith 
3987b3603a34SBarry Smith    Input Parameters:
3988b3603a34SBarry Smith .    ctx - the monitor context
3989b3603a34SBarry Smith 
3990b3603a34SBarry Smith    Level: intermediate
3991b3603a34SBarry Smith 
3992b3603a34SBarry Smith .keywords: TS,  vector, monitor, view
3993b3603a34SBarry Smith 
3994b3603a34SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolutionODE()
3995b3603a34SBarry Smith @*/
3996b3603a34SBarry Smith PetscErrorCode  TSMonitorSolutionODEDestroy(PetscDrawLG *lg)
3997b3603a34SBarry Smith {
3998b3603a34SBarry Smith   PetscDraw      draw;
3999b3603a34SBarry Smith   PetscErrorCode ierr;
4000b3603a34SBarry Smith 
4001b3603a34SBarry Smith   PetscFunctionBegin;
4002b3603a34SBarry Smith   ierr = PetscDrawLGGetDraw(*lg,&draw);CHKERRQ(ierr);
4003b3603a34SBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
4004b3603a34SBarry Smith   ierr = PetscDrawLGDestroy(lg);CHKERRQ(ierr);
4005b3603a34SBarry Smith   PetscFunctionReturn(0);
4006b3603a34SBarry Smith }
4007b3603a34SBarry Smith 
4008b3603a34SBarry Smith #undef __FUNCT__
4009b3603a34SBarry Smith #define __FUNCT__ "TSMonitorSolutionODECreate"
4010b3603a34SBarry Smith /*@C
4011b3603a34SBarry Smith    TSMonitorSolutionODECreate - Creates the monitor context for TSMonitorSolutionODE()
4012b3603a34SBarry Smith 
4013b3603a34SBarry Smith    Collective on TS
4014b3603a34SBarry Smith 
4015b3603a34SBarry Smith    Input Parameter:
4016b3603a34SBarry Smith +    comm - MPI communicator
4017b3603a34SBarry Smith .    N - number of components in the solution vector
4018b3603a34SBarry Smith -
4019b3603a34SBarry Smith 
4020b3603a34SBarry Smith    Output Patameter:
4021b3603a34SBarry Smith .    ctx - the monitor context
4022b3603a34SBarry Smith 
4023b3603a34SBarry Smith    Level: intermediate
4024b3603a34SBarry Smith 
4025b3603a34SBarry Smith .keywords: TS,  vector, monitor, view
4026b3603a34SBarry Smith 
4027b3603a34SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
4028b3603a34SBarry Smith @*/
4029b3603a34SBarry Smith PetscErrorCode  TSMonitorSolutionODECreate(MPI_Comm comm,PetscInt N,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
4030b3603a34SBarry Smith {
4031b3603a34SBarry Smith   PetscDraw      win;
4032b3603a34SBarry Smith   PetscErrorCode ierr;
4033b3603a34SBarry Smith   PetscDrawAxis  axis;
4034b3603a34SBarry Smith 
4035b3603a34SBarry Smith   PetscFunctionBegin;
4036b3603a34SBarry Smith   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr);
4037b3603a34SBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
4038b3603a34SBarry Smith   ierr = PetscDrawLGCreate(win,N,draw);CHKERRQ(ierr);
4039b3603a34SBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
4040b3603a34SBarry Smith   ierr = PetscDrawLGGetAxis(*draw,&axis);CHKERRQ(ierr);
4041b3603a34SBarry Smith   ierr = PetscDrawAxisSetLabels(axis,"Solution","Time","Solution");CHKERRQ(ierr);
4042b3603a34SBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
4043b3603a34SBarry Smith   PetscFunctionReturn(0);
4044b3603a34SBarry Smith }
4045*ef20d060SBarry Smith 
4046*ef20d060SBarry Smith #undef __FUNCT__
4047*ef20d060SBarry Smith #define __FUNCT__ "TSMonitorErrorODE"
4048*ef20d060SBarry Smith /*@C
4049*ef20d060SBarry Smith    TSMonitorErrorODE - Monitors progress of the TS solvers by plotting each component of the solution vector
4050*ef20d060SBarry Smith        in a time based line graph
4051*ef20d060SBarry Smith 
4052*ef20d060SBarry Smith    Collective on TS
4053*ef20d060SBarry Smith 
4054*ef20d060SBarry Smith    Input Parameters:
4055*ef20d060SBarry Smith +  ts - the TS context
4056*ef20d060SBarry Smith .  step - current time-step
4057*ef20d060SBarry Smith .  ptime - current time
4058*ef20d060SBarry Smith -  lg - a line graph object
4059*ef20d060SBarry Smith 
4060*ef20d060SBarry Smith    Level: intermediate
4061*ef20d060SBarry Smith 
4062*ef20d060SBarry Smith     Notes: only for sequential solves
4063*ef20d060SBarry Smith 
4064*ef20d060SBarry Smith .keywords: TS,  vector, monitor, view
4065*ef20d060SBarry Smith 
4066*ef20d060SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
4067*ef20d060SBarry Smith @*/
4068*ef20d060SBarry Smith PetscErrorCode  TSMonitorErrorODE(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
4069*ef20d060SBarry Smith {
4070*ef20d060SBarry Smith   PetscErrorCode    ierr;
4071*ef20d060SBarry Smith   PetscDrawLG       lg = (PetscDrawLG)dummy;
4072*ef20d060SBarry Smith   const PetscScalar *yy;
4073*ef20d060SBarry Smith   Vec               y;
4074*ef20d060SBarry Smith 
4075*ef20d060SBarry Smith   PetscFunctionBegin;
4076*ef20d060SBarry Smith   if (!step) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
4077*ef20d060SBarry Smith   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
4078*ef20d060SBarry Smith   ierr = TSComputeSolutionFunction(ts,ptime,y);CHKERRQ(ierr);
4079*ef20d060SBarry Smith   ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
4080*ef20d060SBarry Smith   ierr = VecGetArrayRead(y,&yy);CHKERRQ(ierr);
4081*ef20d060SBarry Smith   ierr = PetscDrawLGAddCommonPoint(lg,ptime,yy);CHKERRQ(ierr);
4082*ef20d060SBarry Smith   ierr = VecRestoreArrayRead(y,&yy);CHKERRQ(ierr);
4083*ef20d060SBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
4084*ef20d060SBarry Smith   ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
4085*ef20d060SBarry Smith   PetscFunctionReturn(0);
4086*ef20d060SBarry Smith }
4087*ef20d060SBarry Smith 
4088*ef20d060SBarry Smith #undef __FUNCT__
4089*ef20d060SBarry Smith #define __FUNCT__ "TSMonitorErrorODEDestroy"
4090*ef20d060SBarry Smith /*@C
4091*ef20d060SBarry Smith    TSMonitorErrorODEDestroy - Destroys the monitor context for TSMonitorErrorODE()
4092*ef20d060SBarry Smith 
4093*ef20d060SBarry Smith    Collective on TS
4094*ef20d060SBarry Smith 
4095*ef20d060SBarry Smith    Input Parameters:
4096*ef20d060SBarry Smith .    ctx - the monitor context
4097*ef20d060SBarry Smith 
4098*ef20d060SBarry Smith    Level: intermediate
4099*ef20d060SBarry Smith 
4100*ef20d060SBarry Smith .keywords: TS,  vector, monitor, view
4101*ef20d060SBarry Smith 
4102*ef20d060SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorErrorODE()
4103*ef20d060SBarry Smith @*/
4104*ef20d060SBarry Smith PetscErrorCode  TSMonitorErrorODEDestroy(PetscDrawLG *lg)
4105*ef20d060SBarry Smith {
4106*ef20d060SBarry Smith   PetscDraw      draw;
4107*ef20d060SBarry Smith   PetscErrorCode ierr;
4108*ef20d060SBarry Smith 
4109*ef20d060SBarry Smith   PetscFunctionBegin;
4110*ef20d060SBarry Smith   ierr = PetscDrawLGGetDraw(*lg,&draw);CHKERRQ(ierr);
4111*ef20d060SBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
4112*ef20d060SBarry Smith   ierr = PetscDrawLGDestroy(lg);CHKERRQ(ierr);
4113*ef20d060SBarry Smith   PetscFunctionReturn(0);
4114*ef20d060SBarry Smith }
4115*ef20d060SBarry Smith 
4116*ef20d060SBarry Smith #undef __FUNCT__
4117*ef20d060SBarry Smith #define __FUNCT__ "TSMonitorErrorODECreate"
4118*ef20d060SBarry Smith /*@C
4119*ef20d060SBarry Smith    TSMonitorErrorODECreate - Creates the monitor context for TSMonitorErrorODE()
4120*ef20d060SBarry Smith 
4121*ef20d060SBarry Smith    Collective on TS
4122*ef20d060SBarry Smith 
4123*ef20d060SBarry Smith    Input Parameter:
4124*ef20d060SBarry Smith +    comm - MPI communicator
4125*ef20d060SBarry Smith .    N - number of components in the solution vector
4126*ef20d060SBarry Smith -
4127*ef20d060SBarry Smith 
4128*ef20d060SBarry Smith    Output Patameter:
4129*ef20d060SBarry Smith .    ctx - the monitor context
4130*ef20d060SBarry Smith 
4131*ef20d060SBarry Smith    Level: intermediate
4132*ef20d060SBarry Smith 
4133*ef20d060SBarry Smith .keywords: TS,  vector, monitor, view
4134*ef20d060SBarry Smith 
4135*ef20d060SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorError()
4136*ef20d060SBarry Smith @*/
4137*ef20d060SBarry Smith PetscErrorCode  TSMonitorErrorODECreate(MPI_Comm comm,PetscInt N,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
4138*ef20d060SBarry Smith {
4139*ef20d060SBarry Smith   PetscDraw      win;
4140*ef20d060SBarry Smith   PetscErrorCode ierr;
4141*ef20d060SBarry Smith   PetscDrawAxis  axis;
4142*ef20d060SBarry Smith 
4143*ef20d060SBarry Smith   PetscFunctionBegin;
4144*ef20d060SBarry Smith   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr);
4145*ef20d060SBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
4146*ef20d060SBarry Smith   ierr = PetscDrawLGCreate(win,N,draw);CHKERRQ(ierr);
4147*ef20d060SBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
4148*ef20d060SBarry Smith   ierr = PetscDrawLGGetAxis(*draw,&axis);CHKERRQ(ierr);
4149*ef20d060SBarry Smith   ierr = PetscDrawAxisSetLabels(axis,"Error","Time","Error");CHKERRQ(ierr);
4150*ef20d060SBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
4151*ef20d060SBarry Smith   PetscFunctionReturn(0);
4152*ef20d060SBarry Smith }
4153