xref: /petsc/src/ts/interface/ts.c (revision b8123dae7f872ca50579e8005f6d2cd0e4534f13)
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;
80bdad233fSMatthew Knepley 
81bdad233fSMatthew Knepley   PetscFunctionBegin;
820700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
833194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr);
84a43b19c4SJed Brown     /* Handle TS type options */
85a43b19c4SJed Brown     ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr);
86bdad233fSMatthew Knepley 
87bdad233fSMatthew Knepley     /* Handle generic TS options */
88bdad233fSMatthew Knepley     ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr);
893bca7d26SBarry Smith     ierr = PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr);
90d8cd7023SBarry Smith     ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);CHKERRQ(ierr);
91d8cd7023SBarry Smith     ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&ts->time_step,PETSC_NULL);CHKERRQ(ierr);
92a43b19c4SJed Brown     opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time;
93a43b19c4SJed Brown     ierr = PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);CHKERRQ(ierr);
94461e00b3SSean Farley     if (flg) {ierr = TSSetExactFinalTime(ts,opt);CHKERRQ(ierr);}
95cef5090cSJed 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);
96cef5090cSJed 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);
97cef5090cSJed Brown     ierr = PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr);
981c3436cfSJed Brown     ierr = PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,PETSC_NULL);CHKERRQ(ierr);
991c3436cfSJed Brown     ierr = PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,PETSC_NULL);CHKERRQ(ierr);
100bdad233fSMatthew Knepley 
101bdad233fSMatthew Knepley     /* Monitor options */
102a6570f20SBarry Smith     ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
103eabae89aSBarry Smith     if (flg) {
104649052a6SBarry Smith       ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr);
105649052a6SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
106bdad233fSMatthew Knepley     }
1075180491cSLisandro Dalcin     ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1085180491cSLisandro Dalcin     if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);}
1095180491cSLisandro Dalcin 
11090d69ab7SBarry Smith     opt  = PETSC_FALSE;
111acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
112a7cc72afSBarry Smith     if (opt) {
113a6570f20SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
114bdad233fSMatthew Knepley     }
11590d69ab7SBarry Smith     opt  = PETSC_FALSE;
116acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr);
117a7cc72afSBarry Smith     if (opt) {
1186083293cSBarry Smith       void *ctx;
1196083293cSBarry Smith       ierr = TSMonitorSolutionCreate(ts,PETSC_NULL,&ctx);CHKERRQ(ierr);
1206083293cSBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);CHKERRQ(ierr);
121bdad233fSMatthew Knepley     }
122fb1732b5SBarry Smith     opt  = PETSC_FALSE;
123fb1732b5SBarry Smith     ierr = PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
124fb1732b5SBarry Smith     if (flg) {
125fb1732b5SBarry Smith       PetscViewer ctx;
126fb1732b5SBarry Smith       if (monfilename[0]) {
127fb1732b5SBarry Smith         ierr = PetscViewerBinaryOpen(((PetscObject)ts)->comm,monfilename,FILE_MODE_WRITE,&ctx);CHKERRQ(ierr);
128fb1732b5SBarry Smith       } else {
129fb1732b5SBarry Smith         ctx = PETSC_VIEWER_BINARY_(((PetscObject)ts)->comm);
130fb1732b5SBarry Smith       }
131fb1732b5SBarry Smith       ierr = TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
132fb1732b5SBarry Smith     }
133ed81e22dSJed Brown     opt  = PETSC_FALSE;
134ed81e22dSJed 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);
135ed81e22dSJed Brown     if (flg) {
136ed81e22dSJed Brown       const char *ptr,*ptr2;
137ed81e22dSJed Brown       char *filetemplate;
138ed81e22dSJed Brown       if (!monfilename[0]) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
139ed81e22dSJed Brown       /* Do some cursory validation of the input. */
140ed81e22dSJed Brown       ierr = PetscStrstr(monfilename,"%",(char**)&ptr);CHKERRQ(ierr);
141ed81e22dSJed Brown       if (!ptr) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts");
142ed81e22dSJed Brown       for (ptr++ ; ptr && *ptr; ptr++) {
143ed81e22dSJed Brown         ierr = PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);CHKERRQ(ierr);
144ed81e22dSJed 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");
145ed81e22dSJed Brown         if (ptr2) break;
146ed81e22dSJed Brown       }
147ed81e22dSJed Brown       ierr = PetscStrallocpy(monfilename,&filetemplate);CHKERRQ(ierr);
148ed81e22dSJed Brown       ierr = TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);CHKERRQ(ierr);
149ed81e22dSJed Brown     }
150bdad233fSMatthew Knepley 
1511c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
1521c3436cfSJed Brown     ierr = TSAdaptSetFromOptions(adapt);CHKERRQ(ierr);
1531c3436cfSJed Brown 
154d52bd9f3SBarry Smith     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
155d52bd9f3SBarry Smith     if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);}
156d52bd9f3SBarry Smith 
157bdad233fSMatthew Knepley     /* Handle specific TS options */
158abc0a331SBarry Smith     if (ts->ops->setfromoptions) {
159bdad233fSMatthew Knepley       ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr);
160bdad233fSMatthew Knepley     }
1615d973c19SBarry Smith 
1625d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
1635d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr);
164bdad233fSMatthew Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
165bdad233fSMatthew Knepley   PetscFunctionReturn(0);
166bdad233fSMatthew Knepley }
167bdad233fSMatthew Knepley 
168bdad233fSMatthew Knepley #undef __FUNCT__
169cdcf91faSSean Farley #undef __FUNCT__
1704a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSJacobian"
171a7a1495cSBarry Smith /*@
1728c385f81SBarry Smith    TSComputeRHSJacobian - Computes the Jacobian matrix that has been
173a7a1495cSBarry Smith       set with TSSetRHSJacobian().
174a7a1495cSBarry Smith 
175a7a1495cSBarry Smith    Collective on TS and Vec
176a7a1495cSBarry Smith 
177a7a1495cSBarry Smith    Input Parameters:
178316643e7SJed Brown +  ts - the TS context
179a7a1495cSBarry Smith .  t - current timestep
180a7a1495cSBarry Smith -  x - input vector
181a7a1495cSBarry Smith 
182a7a1495cSBarry Smith    Output Parameters:
183a7a1495cSBarry Smith +  A - Jacobian matrix
184a7a1495cSBarry Smith .  B - optional preconditioning matrix
185a7a1495cSBarry Smith -  flag - flag indicating matrix structure
186a7a1495cSBarry Smith 
187a7a1495cSBarry Smith    Notes:
188a7a1495cSBarry Smith    Most users should not need to explicitly call this routine, as it
189a7a1495cSBarry Smith    is used internally within the nonlinear solvers.
190a7a1495cSBarry Smith 
19194b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
192a7a1495cSBarry Smith    flag parameter.
193a7a1495cSBarry Smith 
194a7a1495cSBarry Smith    Level: developer
195a7a1495cSBarry Smith 
196a7a1495cSBarry Smith .keywords: SNES, compute, Jacobian, matrix
197a7a1495cSBarry Smith 
19894b7f48cSBarry Smith .seealso:  TSSetRHSJacobian(), KSPSetOperators()
199a7a1495cSBarry Smith @*/
2007087cfbeSBarry Smith PetscErrorCode  TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg)
201a7a1495cSBarry Smith {
202dfbe8321SBarry Smith   PetscErrorCode ierr;
2030e4ef248SJed Brown   PetscInt Xstate;
204a7a1495cSBarry Smith 
205a7a1495cSBarry Smith   PetscFunctionBegin;
2060700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2070700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
208c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,X,3);
2090e4ef248SJed Brown   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
2100e4ef248SJed Brown   if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) {
2110e4ef248SJed Brown     *flg = ts->rhsjacobian.mstructure;
2120e4ef248SJed Brown     PetscFunctionReturn(0);
213f8ede8e7SJed Brown   }
214d90be118SSean Farley 
2153b5f76d0SSean Farley   if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
216d90be118SSean Farley 
2173b5f76d0SSean Farley   if (ts->userops->rhsjacobian) {
218d5ba7fb7SMatthew Knepley     ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
219a7a1495cSBarry Smith     *flg = DIFFERENT_NONZERO_PATTERN;
220a7a1495cSBarry Smith     PetscStackPush("TS user Jacobian function");
2213b5f76d0SSean Farley     ierr = (*ts->userops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr);
222a7a1495cSBarry Smith     PetscStackPop;
223d5ba7fb7SMatthew Knepley     ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
224a7a1495cSBarry Smith     /* make sure user returned a correct Jacobian and preconditioner */
2250700a824SBarry Smith     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
2260700a824SBarry Smith     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
227ef66eb69SBarry Smith   } else {
228214bc6a2SJed Brown     ierr = MatZeroEntries(*A);CHKERRQ(ierr);
229214bc6a2SJed Brown     if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);}
230214bc6a2SJed Brown     *flg = SAME_NONZERO_PATTERN;
231ef66eb69SBarry Smith   }
2320e4ef248SJed Brown   ts->rhsjacobian.time = t;
2330e4ef248SJed Brown   ts->rhsjacobian.X = X;
2340e4ef248SJed Brown   ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr);
2350e4ef248SJed Brown   ts->rhsjacobian.mstructure = *flg;
236a7a1495cSBarry Smith   PetscFunctionReturn(0);
237a7a1495cSBarry Smith }
238a7a1495cSBarry Smith 
2394a2ae208SSatish Balay #undef __FUNCT__
2404a2ae208SSatish Balay #define __FUNCT__ "TSComputeRHSFunction"
241316643e7SJed Brown /*@
242d763cef2SBarry Smith    TSComputeRHSFunction - Evaluates the right-hand-side function.
243d763cef2SBarry Smith 
244316643e7SJed Brown    Collective on TS and Vec
245316643e7SJed Brown 
246316643e7SJed Brown    Input Parameters:
247316643e7SJed Brown +  ts - the TS context
248316643e7SJed Brown .  t - current time
249316643e7SJed Brown -  x - state vector
250316643e7SJed Brown 
251316643e7SJed Brown    Output Parameter:
252316643e7SJed Brown .  y - right hand side
253316643e7SJed Brown 
254316643e7SJed Brown    Note:
255316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
256316643e7SJed Brown    is used internally within the nonlinear solvers.
257316643e7SJed Brown 
258316643e7SJed Brown    Level: developer
259316643e7SJed Brown 
260316643e7SJed Brown .keywords: TS, compute
261316643e7SJed Brown 
262316643e7SJed Brown .seealso: TSSetRHSFunction(), TSComputeIFunction()
263316643e7SJed Brown @*/
264dfbe8321SBarry Smith PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y)
265d763cef2SBarry Smith {
266dfbe8321SBarry Smith   PetscErrorCode ierr;
267d763cef2SBarry Smith 
268d763cef2SBarry Smith   PetscFunctionBegin;
2690700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2700700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2710700a824SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
272d763cef2SBarry Smith 
2733b5f76d0SSean Farley   if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
274d763cef2SBarry Smith 
275d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
2763b5f76d0SSean Farley   if (ts->userops->rhsfunction) {
277d763cef2SBarry Smith     PetscStackPush("TS user right-hand-side function");
2783b5f76d0SSean Farley     ierr = (*ts->userops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr);
279d763cef2SBarry Smith     PetscStackPop;
280214bc6a2SJed Brown   } else {
281214bc6a2SJed Brown     ierr = VecZeroEntries(y);CHKERRQ(ierr);
282b2cd27e8SJed Brown   }
28344a41b28SSean Farley 
284d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr);
285d763cef2SBarry Smith   PetscFunctionReturn(0);
286d763cef2SBarry Smith }
287d763cef2SBarry Smith 
2884a2ae208SSatish Balay #undef __FUNCT__
289214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSVec_Private"
290214bc6a2SJed Brown static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs)
291214bc6a2SJed Brown {
2922dd45cf8SJed Brown   Vec            F;
293214bc6a2SJed Brown   PetscErrorCode ierr;
294214bc6a2SJed Brown 
295214bc6a2SJed Brown   PetscFunctionBegin;
2960da9a4f0SJed Brown   *Frhs = PETSC_NULL;
2972dd45cf8SJed Brown   ierr = TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
298214bc6a2SJed Brown   if (!ts->Frhs) {
2992dd45cf8SJed Brown     ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr);
300214bc6a2SJed Brown   }
301214bc6a2SJed Brown   *Frhs = ts->Frhs;
302214bc6a2SJed Brown   PetscFunctionReturn(0);
303214bc6a2SJed Brown }
304214bc6a2SJed Brown 
305214bc6a2SJed Brown #undef __FUNCT__
306214bc6a2SJed Brown #define __FUNCT__ "TSGetRHSMats_Private"
307214bc6a2SJed Brown static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs)
308214bc6a2SJed Brown {
309214bc6a2SJed Brown   Mat            A,B;
3102dd45cf8SJed Brown   PetscErrorCode ierr;
311214bc6a2SJed Brown 
312214bc6a2SJed Brown   PetscFunctionBegin;
313214bc6a2SJed Brown   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
314214bc6a2SJed Brown   if (Arhs) {
315214bc6a2SJed Brown     if (!ts->Arhs) {
316214bc6a2SJed Brown       ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr);
317214bc6a2SJed Brown     }
318214bc6a2SJed Brown     *Arhs = ts->Arhs;
319214bc6a2SJed Brown   }
320214bc6a2SJed Brown   if (Brhs) {
321214bc6a2SJed Brown     if (!ts->Brhs) {
322214bc6a2SJed Brown       ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr);
323214bc6a2SJed Brown     }
324214bc6a2SJed Brown     *Brhs = ts->Brhs;
325214bc6a2SJed Brown   }
326214bc6a2SJed Brown   PetscFunctionReturn(0);
327214bc6a2SJed Brown }
328214bc6a2SJed Brown 
329214bc6a2SJed Brown #undef __FUNCT__
330316643e7SJed Brown #define __FUNCT__ "TSComputeIFunction"
331316643e7SJed Brown /*@
332316643e7SJed Brown    TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0
333316643e7SJed Brown 
334316643e7SJed Brown    Collective on TS and Vec
335316643e7SJed Brown 
336316643e7SJed Brown    Input Parameters:
337316643e7SJed Brown +  ts - the TS context
338316643e7SJed Brown .  t - current time
339316643e7SJed Brown .  X - state vector
340214bc6a2SJed Brown .  Xdot - time derivative of state vector
341214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate
342316643e7SJed Brown 
343316643e7SJed Brown    Output Parameter:
344316643e7SJed Brown .  Y - right hand side
345316643e7SJed Brown 
346316643e7SJed Brown    Note:
347316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
348316643e7SJed Brown    is used internally within the nonlinear solvers.
349316643e7SJed Brown 
350316643e7SJed Brown    If the user did did not write their equations in implicit form, this
351316643e7SJed Brown    function recasts them in implicit form.
352316643e7SJed Brown 
353316643e7SJed Brown    Level: developer
354316643e7SJed Brown 
355316643e7SJed Brown .keywords: TS, compute
356316643e7SJed Brown 
357316643e7SJed Brown .seealso: TSSetIFunction(), TSComputeRHSFunction()
358316643e7SJed Brown @*/
359214bc6a2SJed Brown PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex)
360316643e7SJed Brown {
361316643e7SJed Brown   PetscErrorCode ierr;
362316643e7SJed Brown 
363316643e7SJed Brown   PetscFunctionBegin;
3640700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3650700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
3660700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
3670700a824SBarry Smith   PetscValidHeaderSpecific(Y,VEC_CLASSID,5);
368316643e7SJed Brown 
369d90be118SSean Farley   if (!ts->userops->rhsfunction && !ts->userops->ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()");
370d90be118SSean Farley 
371316643e7SJed Brown   ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
3723b5f76d0SSean Farley   if (ts->userops->ifunction) {
373316643e7SJed Brown     PetscStackPush("TS user implicit function");
3743b5f76d0SSean Farley     ierr = (*ts->userops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr);
375316643e7SJed Brown     PetscStackPop;
376214bc6a2SJed Brown   }
377214bc6a2SJed Brown   if (imex) {
3783b5f76d0SSean Farley     if (!ts->userops->ifunction) {
3792dd45cf8SJed Brown       ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);
3802dd45cf8SJed Brown     }
3812dd45cf8SJed Brown   } else if (ts->userops->rhsfunction) {
3822dd45cf8SJed Brown     if (ts->userops->ifunction) {
383214bc6a2SJed Brown       Vec Frhs;
384214bc6a2SJed Brown       ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr);
385214bc6a2SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr);
386214bc6a2SJed Brown       ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr);
3872dd45cf8SJed Brown     } else {
3882dd45cf8SJed Brown       ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr);
3892dd45cf8SJed Brown       ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr);
390316643e7SJed Brown     }
3914a6899ffSJed Brown   }
392316643e7SJed Brown   ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr);
393316643e7SJed Brown   PetscFunctionReturn(0);
394316643e7SJed Brown }
395316643e7SJed Brown 
396316643e7SJed Brown #undef __FUNCT__
397316643e7SJed Brown #define __FUNCT__ "TSComputeIJacobian"
398316643e7SJed Brown /*@
399316643e7SJed Brown    TSComputeIJacobian - Evaluates the Jacobian of the DAE
400316643e7SJed Brown 
401316643e7SJed Brown    Collective on TS and Vec
402316643e7SJed Brown 
403316643e7SJed Brown    Input
404316643e7SJed Brown       Input Parameters:
405316643e7SJed Brown +  ts - the TS context
406316643e7SJed Brown .  t - current timestep
407316643e7SJed Brown .  X - state vector
408316643e7SJed Brown .  Xdot - time derivative of state vector
409214bc6a2SJed Brown .  shift - shift to apply, see note below
410214bc6a2SJed Brown -  imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
411316643e7SJed Brown 
412316643e7SJed Brown    Output Parameters:
413316643e7SJed Brown +  A - Jacobian matrix
414316643e7SJed Brown .  B - optional preconditioning matrix
415316643e7SJed Brown -  flag - flag indicating matrix structure
416316643e7SJed Brown 
417316643e7SJed Brown    Notes:
418316643e7SJed Brown    If F(t,X,Xdot)=0 is the DAE, the required Jacobian is
419316643e7SJed Brown 
420316643e7SJed Brown    dF/dX + shift*dF/dXdot
421316643e7SJed Brown 
422316643e7SJed Brown    Most users should not need to explicitly call this routine, as it
423316643e7SJed Brown    is used internally within the nonlinear solvers.
424316643e7SJed Brown 
425316643e7SJed Brown    Level: developer
426316643e7SJed Brown 
427316643e7SJed Brown .keywords: TS, compute, Jacobian, matrix
428316643e7SJed Brown 
429316643e7SJed Brown .seealso:  TSSetIJacobian()
43063495f91SJed Brown @*/
431214bc6a2SJed Brown PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex)
432316643e7SJed Brown {
4330026cea9SSean Farley   PetscInt Xstate, Xdotstate;
434316643e7SJed Brown   PetscErrorCode ierr;
435316643e7SJed Brown 
436316643e7SJed Brown   PetscFunctionBegin;
4370700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4380700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
4390700a824SBarry Smith   PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4);
440316643e7SJed Brown   PetscValidPointer(A,6);
4410700a824SBarry Smith   PetscValidHeaderSpecific(*A,MAT_CLASSID,6);
442316643e7SJed Brown   PetscValidPointer(B,7);
4430700a824SBarry Smith   PetscValidHeaderSpecific(*B,MAT_CLASSID,7);
444316643e7SJed Brown   PetscValidPointer(flg,8);
4450026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr);
4460026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr);
4470026cea9SSean 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))) {
4480026cea9SSean Farley     *flg = ts->ijacobian.mstructure;
4490026cea9SSean Farley     ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr);
4500026cea9SSean Farley     PetscFunctionReturn(0);
4510026cea9SSean Farley   }
452316643e7SJed Brown 
4533b5f76d0SSean Farley   if (!ts->userops->rhsjacobian && !ts->userops->ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()");
454316643e7SJed Brown 
4554e684422SJed Brown   *flg = SAME_NONZERO_PATTERN;  /* In case we're solving a linear problem in which case it wouldn't get initialized below. */
456316643e7SJed Brown   ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
4573b5f76d0SSean Farley   if (ts->userops->ijacobian) {
4582dd45cf8SJed Brown     *flg = DIFFERENT_NONZERO_PATTERN;
459316643e7SJed Brown     PetscStackPush("TS user implicit Jacobian");
4603b5f76d0SSean Farley     ierr = (*ts->userops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr);
461316643e7SJed Brown     PetscStackPop;
462214bc6a2SJed Brown     /* make sure user returned a correct Jacobian and preconditioner */
463214bc6a2SJed Brown     PetscValidHeaderSpecific(*A,MAT_CLASSID,4);
464214bc6a2SJed Brown     PetscValidHeaderSpecific(*B,MAT_CLASSID,5);
4654a6899ffSJed Brown   }
466214bc6a2SJed Brown   if (imex) {
4673b5f76d0SSean Farley     if (!ts->userops->ijacobian) {  /* system was written as Xdot = F(t,X) */
468214bc6a2SJed Brown       ierr = MatZeroEntries(*A);CHKERRQ(ierr);
4692dd45cf8SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
470214bc6a2SJed Brown       if (*A != *B) {
471214bc6a2SJed Brown         ierr = MatZeroEntries(*B);CHKERRQ(ierr);
472214bc6a2SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
473214bc6a2SJed Brown       }
474214bc6a2SJed Brown       *flg = SAME_PRECONDITIONER;
475214bc6a2SJed Brown     }
476214bc6a2SJed Brown   } else {
4773b5f76d0SSean Farley     if (!ts->userops->ijacobian) {
478214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr);
479214bc6a2SJed Brown       ierr = MatScale(*A,-1);CHKERRQ(ierr);
480214bc6a2SJed Brown       ierr = MatShift(*A,shift);CHKERRQ(ierr);
481316643e7SJed Brown       if (*A != *B) {
482316643e7SJed Brown         ierr = MatScale(*B,-1);CHKERRQ(ierr);
483316643e7SJed Brown         ierr = MatShift(*B,shift);CHKERRQ(ierr);
484316643e7SJed Brown       }
4853b5f76d0SSean Farley     } else if (ts->userops->rhsjacobian) {
486214bc6a2SJed Brown       Mat Arhs,Brhs;
487214bc6a2SJed Brown       MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN;
488214bc6a2SJed Brown       ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
489214bc6a2SJed Brown       ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
490214bc6a2SJed Brown       axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN;
491214bc6a2SJed Brown       ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr);
492214bc6a2SJed Brown       if (*A != *B) {
493214bc6a2SJed Brown         ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr);
494214bc6a2SJed Brown       }
495214bc6a2SJed Brown       *flg = PetscMin(*flg,flg2);
496214bc6a2SJed Brown     }
497316643e7SJed Brown   }
4980026cea9SSean Farley 
4990026cea9SSean Farley   ts->ijacobian.time = t;
5000026cea9SSean Farley   ts->ijacobian.X = X;
5010026cea9SSean Farley   ts->ijacobian.Xdot = Xdot;
5020026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr);
5030026cea9SSean Farley   ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr);
5040026cea9SSean Farley   ts->ijacobian.shift = shift;
5050026cea9SSean Farley   ts->ijacobian.imex = imex;
5060026cea9SSean Farley   ts->ijacobian.mstructure = *flg;
507316643e7SJed Brown   ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr);
508316643e7SJed Brown   PetscFunctionReturn(0);
509316643e7SJed Brown }
510316643e7SJed Brown 
511316643e7SJed Brown #undef __FUNCT__
5124a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSFunction"
513d763cef2SBarry Smith /*@C
514d763cef2SBarry Smith     TSSetRHSFunction - Sets the routine for evaluating the function,
515d763cef2SBarry Smith     F(t,u), where U_t = F(t,u).
516d763cef2SBarry Smith 
5173f9fe445SBarry Smith     Logically Collective on TS
518d763cef2SBarry Smith 
519d763cef2SBarry Smith     Input Parameters:
520d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
521ca94891dSJed Brown .   r - vector to put the computed right hand side (or PETSC_NULL to have it created)
522d763cef2SBarry Smith .   f - routine for evaluating the right-hand-side function
523d763cef2SBarry Smith -   ctx - [optional] user-defined context for private data for the
524d763cef2SBarry Smith           function evaluation routine (may be PETSC_NULL)
525d763cef2SBarry Smith 
526d763cef2SBarry Smith     Calling sequence of func:
52787828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Vec F,void *ctx);
528d763cef2SBarry Smith 
529d763cef2SBarry Smith +   t - current timestep
530d763cef2SBarry Smith .   u - input vector
531d763cef2SBarry Smith .   F - function vector
532d763cef2SBarry Smith -   ctx - [optional] user-defined function context
533d763cef2SBarry Smith 
534d763cef2SBarry Smith     Level: beginner
535d763cef2SBarry Smith 
536d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, function
537d763cef2SBarry Smith 
538d6cbdb99SBarry Smith .seealso: TSSetRHSJacobian(), TSSetIJacobian()
539d763cef2SBarry Smith @*/
540089b2837SJed Brown PetscErrorCode  TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx)
541d763cef2SBarry Smith {
542089b2837SJed Brown   PetscErrorCode ierr;
543089b2837SJed Brown   SNES           snes;
544e856ceecSJed Brown   Vec            ralloc = PETSC_NULL;
545d763cef2SBarry Smith 
546089b2837SJed Brown   PetscFunctionBegin;
5470700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
548ca94891dSJed Brown   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
5493b5f76d0SSean Farley   if (f)   ts->userops->rhsfunction = f;
5504e684422SJed Brown   if (ctx) ts->funP                 = ctx;
551089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
552e856ceecSJed Brown   if (!r && !ts->dm && ts->vec_sol) {
553e856ceecSJed Brown     ierr = VecDuplicate(ts->vec_sol,&ralloc);CHKERRQ(ierr);
554e856ceecSJed Brown     r = ralloc;
555e856ceecSJed Brown   }
556089b2837SJed Brown   ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr);
557e856ceecSJed Brown   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
558d763cef2SBarry Smith   PetscFunctionReturn(0);
559d763cef2SBarry Smith }
560d763cef2SBarry Smith 
5614a2ae208SSatish Balay #undef __FUNCT__
5624a2ae208SSatish Balay #define __FUNCT__ "TSSetRHSJacobian"
563d763cef2SBarry Smith /*@C
564d763cef2SBarry Smith    TSSetRHSJacobian - Sets the function to compute the Jacobian of F,
565d763cef2SBarry Smith    where U_t = F(U,t), as well as the location to store the matrix.
566d763cef2SBarry Smith 
5673f9fe445SBarry Smith    Logically Collective on TS
568d763cef2SBarry Smith 
569d763cef2SBarry Smith    Input Parameters:
570d763cef2SBarry Smith +  ts  - the TS context obtained from TSCreate()
571d763cef2SBarry Smith .  A   - Jacobian matrix
572d763cef2SBarry Smith .  B   - preconditioner matrix (usually same as A)
573d763cef2SBarry Smith .  f   - the Jacobian evaluation routine
574d763cef2SBarry Smith -  ctx - [optional] user-defined context for private data for the
575d763cef2SBarry Smith          Jacobian evaluation routine (may be PETSC_NULL)
576d763cef2SBarry Smith 
577d763cef2SBarry Smith    Calling sequence of func:
57887828ca2SBarry Smith $     func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx);
579d763cef2SBarry Smith 
580d763cef2SBarry Smith +  t - current timestep
581d763cef2SBarry Smith .  u - input vector
582d763cef2SBarry Smith .  A - matrix A, where U_t = A(t)u
583d763cef2SBarry Smith .  B - preconditioner matrix, usually the same as A
584d763cef2SBarry Smith .  flag - flag indicating information about the preconditioner matrix
58594b7f48cSBarry Smith           structure (same as flag in KSPSetOperators())
586d763cef2SBarry Smith -  ctx - [optional] user-defined context for matrix evaluation routine
587d763cef2SBarry Smith 
588d763cef2SBarry Smith    Notes:
58994b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
590d763cef2SBarry Smith    output parameter in the routine func().  Be sure to read this information!
591d763cef2SBarry Smith 
592d763cef2SBarry Smith    The routine func() takes Mat * as the matrix arguments rather than Mat.
593d763cef2SBarry Smith    This allows the matrix evaluation routine to replace A and/or B with a
59456335db2SHong Zhang    completely new matrix structure (not just different matrix elements)
595d763cef2SBarry Smith    when appropriate, for instance, if the nonzero structure is changing
596d763cef2SBarry Smith    throughout the global iterations.
597d763cef2SBarry Smith 
598d763cef2SBarry Smith    Level: beginner
599d763cef2SBarry Smith 
600d763cef2SBarry Smith .keywords: TS, timestep, set, right-hand-side, Jacobian
601d763cef2SBarry Smith 
602ab637aeaSJed Brown .seealso: SNESDefaultComputeJacobianColor(), TSSetRHSFunction()
603d763cef2SBarry Smith 
604d763cef2SBarry Smith @*/
605089b2837SJed Brown PetscErrorCode  TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx)
606d763cef2SBarry Smith {
607277b19d0SLisandro Dalcin   PetscErrorCode ierr;
608089b2837SJed Brown   SNES           snes;
609277b19d0SLisandro Dalcin 
610d763cef2SBarry Smith   PetscFunctionBegin;
6110700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
612277b19d0SLisandro Dalcin   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
613277b19d0SLisandro Dalcin   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
614277b19d0SLisandro Dalcin   if (A) PetscCheckSameComm(ts,1,A,2);
615277b19d0SLisandro Dalcin   if (B) PetscCheckSameComm(ts,1,B,3);
616d763cef2SBarry Smith 
6173b5f76d0SSean Farley   if (f)   ts->userops->rhsjacobian = f;
618277b19d0SLisandro Dalcin   if (ctx) ts->jacP                 = ctx;
619089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
620dabe61deSJed Brown   if (!ts->userops->ijacobian) {
621089b2837SJed Brown     ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
6220e4ef248SJed Brown   }
6230e4ef248SJed Brown   if (A) {
6240e4ef248SJed Brown     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
6250e4ef248SJed Brown     ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
6260e4ef248SJed Brown     ts->Arhs = A;
6270e4ef248SJed Brown   }
6280e4ef248SJed Brown   if (B) {
6290e4ef248SJed Brown     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
6300e4ef248SJed Brown     ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
6310e4ef248SJed Brown     ts->Brhs = B;
6320e4ef248SJed Brown   }
633d763cef2SBarry Smith   PetscFunctionReturn(0);
634d763cef2SBarry Smith }
635d763cef2SBarry Smith 
636316643e7SJed Brown 
637316643e7SJed Brown #undef __FUNCT__
638316643e7SJed Brown #define __FUNCT__ "TSSetIFunction"
639316643e7SJed Brown /*@C
640316643e7SJed Brown    TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved.
641316643e7SJed Brown 
6423f9fe445SBarry Smith    Logically Collective on TS
643316643e7SJed Brown 
644316643e7SJed Brown    Input Parameters:
645316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
646ca94891dSJed Brown .  r   - vector to hold the residual (or PETSC_NULL to have it created internally)
647316643e7SJed Brown .  f   - the function evaluation routine
648316643e7SJed Brown -  ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL)
649316643e7SJed Brown 
650316643e7SJed Brown    Calling sequence of f:
651316643e7SJed Brown $  f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx);
652316643e7SJed Brown 
653316643e7SJed Brown +  t   - time at step/stage being solved
654316643e7SJed Brown .  u   - state vector
655316643e7SJed Brown .  u_t - time derivative of state vector
656316643e7SJed Brown .  F   - function vector
657316643e7SJed Brown -  ctx - [optional] user-defined context for matrix evaluation routine
658316643e7SJed Brown 
659316643e7SJed Brown    Important:
660d6cbdb99SBarry Smith    The user MUST call either this routine, TSSetRHSFunction().  This routine must be used when not solving an ODE, for example a DAE.
661316643e7SJed Brown 
662316643e7SJed Brown    Level: beginner
663316643e7SJed Brown 
664316643e7SJed Brown .keywords: TS, timestep, set, DAE, Jacobian
665316643e7SJed Brown 
666d6cbdb99SBarry Smith .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian()
667316643e7SJed Brown @*/
668089b2837SJed Brown PetscErrorCode  TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx)
669316643e7SJed Brown {
670089b2837SJed Brown   PetscErrorCode ierr;
671089b2837SJed Brown   SNES           snes;
672e856ceecSJed Brown   Vec            resalloc = PETSC_NULL;
673316643e7SJed Brown 
674316643e7SJed Brown   PetscFunctionBegin;
6750700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
676ca94891dSJed Brown   if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2);
6773b5f76d0SSean Farley   if (f)   ts->userops->ifunction = f;
678089b2837SJed Brown   if (ctx) ts->funP           = ctx;
679089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
680e856ceecSJed Brown   if (!res && !ts->dm && ts->vec_sol) {
681e856ceecSJed Brown     ierr = VecDuplicate(ts->vec_sol,&resalloc);CHKERRQ(ierr);
682e856ceecSJed Brown     res = resalloc;
683e856ceecSJed Brown   }
684089b2837SJed Brown   ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr);
685e856ceecSJed Brown   ierr = VecDestroy(&resalloc);CHKERRQ(ierr);
686089b2837SJed Brown   PetscFunctionReturn(0);
687089b2837SJed Brown }
688089b2837SJed Brown 
689089b2837SJed Brown #undef __FUNCT__
690089b2837SJed Brown #define __FUNCT__ "TSGetIFunction"
691089b2837SJed Brown /*@C
692089b2837SJed Brown    TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it.
693089b2837SJed Brown 
694089b2837SJed Brown    Not Collective
695089b2837SJed Brown 
696089b2837SJed Brown    Input Parameter:
697089b2837SJed Brown .  ts - the TS context
698089b2837SJed Brown 
699089b2837SJed Brown    Output Parameter:
700089b2837SJed Brown +  r - vector to hold residual (or PETSC_NULL)
701089b2837SJed Brown .  func - the function to compute residual (or PETSC_NULL)
702089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
703089b2837SJed Brown 
704089b2837SJed Brown    Level: advanced
705089b2837SJed Brown 
706089b2837SJed Brown .keywords: TS, nonlinear, get, function
707089b2837SJed Brown 
708089b2837SJed Brown .seealso: TSSetIFunction(), SNESGetFunction()
709089b2837SJed Brown @*/
710089b2837SJed Brown PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx)
711089b2837SJed Brown {
712089b2837SJed Brown   PetscErrorCode ierr;
713089b2837SJed Brown   SNES snes;
714089b2837SJed Brown 
715089b2837SJed Brown   PetscFunctionBegin;
716089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
717089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
718089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
7193b5f76d0SSean Farley   if (func) *func = ts->userops->ifunction;
720089b2837SJed Brown   if (ctx)  *ctx  = ts->funP;
721089b2837SJed Brown   PetscFunctionReturn(0);
722089b2837SJed Brown }
723089b2837SJed Brown 
724089b2837SJed Brown #undef __FUNCT__
725089b2837SJed Brown #define __FUNCT__ "TSGetRHSFunction"
726089b2837SJed Brown /*@C
727089b2837SJed Brown    TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it.
728089b2837SJed Brown 
729089b2837SJed Brown    Not Collective
730089b2837SJed Brown 
731089b2837SJed Brown    Input Parameter:
732089b2837SJed Brown .  ts - the TS context
733089b2837SJed Brown 
734089b2837SJed Brown    Output Parameter:
735089b2837SJed Brown +  r - vector to hold computed right hand side (or PETSC_NULL)
736089b2837SJed Brown .  func - the function to compute right hand side (or PETSC_NULL)
737089b2837SJed Brown -  ctx - the function context (or PETSC_NULL)
738089b2837SJed Brown 
739089b2837SJed Brown    Level: advanced
740089b2837SJed Brown 
741089b2837SJed Brown .keywords: TS, nonlinear, get, function
742089b2837SJed Brown 
743089b2837SJed Brown .seealso: TSSetRhsfunction(), SNESGetFunction()
744089b2837SJed Brown @*/
745089b2837SJed Brown PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx)
746089b2837SJed Brown {
747089b2837SJed Brown   PetscErrorCode ierr;
748089b2837SJed Brown   SNES snes;
749089b2837SJed Brown 
750089b2837SJed Brown   PetscFunctionBegin;
751089b2837SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
752089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
753089b2837SJed Brown   ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
7543b5f76d0SSean Farley   if (func) *func = ts->userops->rhsfunction;
755089b2837SJed Brown   if (ctx)  *ctx  = ts->funP;
756316643e7SJed Brown   PetscFunctionReturn(0);
757316643e7SJed Brown }
758316643e7SJed Brown 
759316643e7SJed Brown #undef __FUNCT__
760316643e7SJed Brown #define __FUNCT__ "TSSetIJacobian"
761316643e7SJed Brown /*@C
762a4f0a591SBarry Smith    TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function
763a4f0a591SBarry Smith         you provided with TSSetIFunction().
764316643e7SJed Brown 
7653f9fe445SBarry Smith    Logically Collective on TS
766316643e7SJed Brown 
767316643e7SJed Brown    Input Parameters:
768316643e7SJed Brown +  ts  - the TS context obtained from TSCreate()
769316643e7SJed Brown .  A   - Jacobian matrix
770316643e7SJed Brown .  B   - preconditioning matrix for A (may be same as A)
771316643e7SJed Brown .  f   - the Jacobian evaluation routine
772316643e7SJed Brown -  ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL)
773316643e7SJed Brown 
774316643e7SJed Brown    Calling sequence of f:
7751b4a444bSJed Brown $  f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx);
776316643e7SJed Brown 
777316643e7SJed Brown +  t    - time at step/stage being solved
7781b4a444bSJed Brown .  U    - state vector
7791b4a444bSJed Brown .  U_t  - time derivative of state vector
780316643e7SJed Brown .  a    - shift
7811b4a444bSJed Brown .  A    - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t
782316643e7SJed Brown .  B    - preconditioning matrix for A, may be same as A
783316643e7SJed Brown .  flag - flag indicating information about the preconditioner matrix
784316643e7SJed Brown           structure (same as flag in KSPSetOperators())
785316643e7SJed Brown -  ctx  - [optional] user-defined context for matrix evaluation routine
786316643e7SJed Brown 
787316643e7SJed Brown    Notes:
788316643e7SJed Brown    The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve.
789316643e7SJed Brown 
790a4f0a591SBarry Smith    The matrix dF/dU + a*dF/dU_t you provide turns out to be
791a4f0a591SBarry 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.
792a4f0a591SBarry Smith    The time integrator internally approximates U_t by W+a*U where the positive "shift"
793a4f0a591SBarry Smith    a and vector W depend on the integration method, step size, and past states. For example with
794a4f0a591SBarry Smith    the backward Euler method a = 1/dt and W = -a*U(previous timestep) so
795a4f0a591SBarry Smith    W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt
796a4f0a591SBarry Smith 
797316643e7SJed Brown    Level: beginner
798316643e7SJed Brown 
799316643e7SJed Brown .keywords: TS, timestep, DAE, Jacobian
800316643e7SJed Brown 
801ab637aeaSJed Brown .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESDefaultComputeJacobianColor(), SNESDefaultComputeJacobian()
802316643e7SJed Brown 
803316643e7SJed Brown @*/
8047087cfbeSBarry Smith PetscErrorCode  TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx)
805316643e7SJed Brown {
806316643e7SJed Brown   PetscErrorCode ierr;
807089b2837SJed Brown   SNES           snes;
808316643e7SJed Brown 
809316643e7SJed Brown   PetscFunctionBegin;
8100700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8110700a824SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
8120700a824SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
813316643e7SJed Brown   if (A) PetscCheckSameComm(ts,1,A,2);
814316643e7SJed Brown   if (B) PetscCheckSameComm(ts,1,B,3);
8153b5f76d0SSean Farley   if (f)   ts->userops->ijacobian = f;
816316643e7SJed Brown   if (ctx) ts->jacP           = ctx;
817089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
818089b2837SJed Brown   ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr);
819316643e7SJed Brown   PetscFunctionReturn(0);
820316643e7SJed Brown }
821316643e7SJed Brown 
8224a2ae208SSatish Balay #undef __FUNCT__
8234a2ae208SSatish Balay #define __FUNCT__ "TSView"
8247e2c5f70SBarry Smith /*@C
825d763cef2SBarry Smith     TSView - Prints the TS data structure.
826d763cef2SBarry Smith 
8274c49b128SBarry Smith     Collective on TS
828d763cef2SBarry Smith 
829d763cef2SBarry Smith     Input Parameters:
830d763cef2SBarry Smith +   ts - the TS context obtained from TSCreate()
831d763cef2SBarry Smith -   viewer - visualization context
832d763cef2SBarry Smith 
833d763cef2SBarry Smith     Options Database Key:
834d763cef2SBarry Smith .   -ts_view - calls TSView() at end of TSStep()
835d763cef2SBarry Smith 
836d763cef2SBarry Smith     Notes:
837d763cef2SBarry Smith     The available visualization contexts include
838b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
839b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
840d763cef2SBarry Smith          output where only the first processor opens
841d763cef2SBarry Smith          the file.  All other processors send their
842d763cef2SBarry Smith          data to the first processor to print.
843d763cef2SBarry Smith 
844d763cef2SBarry Smith     The user can open an alternative visualization context with
845b0a32e0cSBarry Smith     PetscViewerASCIIOpen() - output to a specified file.
846d763cef2SBarry Smith 
847d763cef2SBarry Smith     Level: beginner
848d763cef2SBarry Smith 
849d763cef2SBarry Smith .keywords: TS, timestep, view
850d763cef2SBarry Smith 
851b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
852d763cef2SBarry Smith @*/
8537087cfbeSBarry Smith PetscErrorCode  TSView(TS ts,PetscViewer viewer)
854d763cef2SBarry Smith {
855dfbe8321SBarry Smith   PetscErrorCode ierr;
856a313700dSBarry Smith   const TSType   type;
857089b2837SJed Brown   PetscBool      iascii,isstring,isundials;
858d763cef2SBarry Smith 
859d763cef2SBarry Smith   PetscFunctionBegin;
8600700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8613050cee2SBarry Smith   if (!viewer) {
8627adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr);
8633050cee2SBarry Smith   }
8640700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
865c9780b6fSBarry Smith   PetscCheckSameComm(ts,1,viewer,2);
866fd16b177SBarry Smith 
867251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
868251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
86932077d6dSBarry Smith   if (iascii) {
870317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr);
87177431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr);
872a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum time=%G\n",ts->max_time);CHKERRQ(ierr);
873d763cef2SBarry Smith     if (ts->problem_type == TS_NONLINEAR) {
8745ef26d82SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solver iterations=%D\n",ts->snes_its);CHKERRQ(ierr);
875c610991cSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr);
876d763cef2SBarry Smith     }
8775ef26d82SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",ts->ksp_its);CHKERRQ(ierr);
878193ac0bcSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr);
879d52bd9f3SBarry Smith     if (ts->ops->view) {
880d52bd9f3SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
881d52bd9f3SBarry Smith       ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr);
882d52bd9f3SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
883d52bd9f3SBarry Smith     }
8840f5bd95cSBarry Smith   } else if (isstring) {
885a313700dSBarry Smith     ierr = TSGetType(ts,&type);CHKERRQ(ierr);
886b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr);
887d763cef2SBarry Smith   }
888b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
889251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr);
890b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
891d763cef2SBarry Smith   PetscFunctionReturn(0);
892d763cef2SBarry Smith }
893d763cef2SBarry Smith 
894d763cef2SBarry Smith 
8954a2ae208SSatish Balay #undef __FUNCT__
8964a2ae208SSatish Balay #define __FUNCT__ "TSSetApplicationContext"
897b07ff414SBarry Smith /*@
898d763cef2SBarry Smith    TSSetApplicationContext - Sets an optional user-defined context for
899d763cef2SBarry Smith    the timesteppers.
900d763cef2SBarry Smith 
9013f9fe445SBarry Smith    Logically Collective on TS
902d763cef2SBarry Smith 
903d763cef2SBarry Smith    Input Parameters:
904d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
905d763cef2SBarry Smith -  usrP - optional user context
906d763cef2SBarry Smith 
907d763cef2SBarry Smith    Level: intermediate
908d763cef2SBarry Smith 
909d763cef2SBarry Smith .keywords: TS, timestep, set, application, context
910d763cef2SBarry Smith 
911d763cef2SBarry Smith .seealso: TSGetApplicationContext()
912d763cef2SBarry Smith @*/
9137087cfbeSBarry Smith PetscErrorCode  TSSetApplicationContext(TS ts,void *usrP)
914d763cef2SBarry Smith {
915d763cef2SBarry Smith   PetscFunctionBegin;
9160700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
917d763cef2SBarry Smith   ts->user = usrP;
918d763cef2SBarry Smith   PetscFunctionReturn(0);
919d763cef2SBarry Smith }
920d763cef2SBarry Smith 
9214a2ae208SSatish Balay #undef __FUNCT__
9224a2ae208SSatish Balay #define __FUNCT__ "TSGetApplicationContext"
923b07ff414SBarry Smith /*@
924d763cef2SBarry Smith     TSGetApplicationContext - Gets the user-defined context for the
925d763cef2SBarry Smith     timestepper.
926d763cef2SBarry Smith 
927d763cef2SBarry Smith     Not Collective
928d763cef2SBarry Smith 
929d763cef2SBarry Smith     Input Parameter:
930d763cef2SBarry Smith .   ts - the TS context obtained from TSCreate()
931d763cef2SBarry Smith 
932d763cef2SBarry Smith     Output Parameter:
933d763cef2SBarry Smith .   usrP - user context
934d763cef2SBarry Smith 
935d763cef2SBarry Smith     Level: intermediate
936d763cef2SBarry Smith 
937d763cef2SBarry Smith .keywords: TS, timestep, get, application, context
938d763cef2SBarry Smith 
939d763cef2SBarry Smith .seealso: TSSetApplicationContext()
940d763cef2SBarry Smith @*/
941e71120c6SJed Brown PetscErrorCode  TSGetApplicationContext(TS ts,void *usrP)
942d763cef2SBarry Smith {
943d763cef2SBarry Smith   PetscFunctionBegin;
9440700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
945e71120c6SJed Brown   *(void**)usrP = ts->user;
946d763cef2SBarry Smith   PetscFunctionReturn(0);
947d763cef2SBarry Smith }
948d763cef2SBarry Smith 
9494a2ae208SSatish Balay #undef __FUNCT__
9504a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStepNumber"
951d763cef2SBarry Smith /*@
952*b8123daeSJed Brown    TSGetTimeStepNumber - Gets the number of time steps completed.
953d763cef2SBarry Smith 
954d763cef2SBarry Smith    Not Collective
955d763cef2SBarry Smith 
956d763cef2SBarry Smith    Input Parameter:
957d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
958d763cef2SBarry Smith 
959d763cef2SBarry Smith    Output Parameter:
960*b8123daeSJed Brown .  iter - number of steps completed so far
961d763cef2SBarry Smith 
962d763cef2SBarry Smith    Level: intermediate
963d763cef2SBarry Smith 
964d763cef2SBarry Smith .keywords: TS, timestep, get, iteration, number
965*b8123daeSJed Brown .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStep()
966d763cef2SBarry Smith @*/
9677087cfbeSBarry Smith PetscErrorCode  TSGetTimeStepNumber(TS ts,PetscInt* iter)
968d763cef2SBarry Smith {
969d763cef2SBarry Smith   PetscFunctionBegin;
9700700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
9714482741eSBarry Smith   PetscValidIntPointer(iter,2);
972d763cef2SBarry Smith   *iter = ts->steps;
973d763cef2SBarry Smith   PetscFunctionReturn(0);
974d763cef2SBarry Smith }
975d763cef2SBarry Smith 
9764a2ae208SSatish Balay #undef __FUNCT__
9774a2ae208SSatish Balay #define __FUNCT__ "TSSetInitialTimeStep"
978d763cef2SBarry Smith /*@
979d763cef2SBarry Smith    TSSetInitialTimeStep - Sets the initial timestep to be used,
980d763cef2SBarry Smith    as well as the initial time.
981d763cef2SBarry Smith 
9823f9fe445SBarry Smith    Logically Collective on TS
983d763cef2SBarry Smith 
984d763cef2SBarry Smith    Input Parameters:
985d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
986d763cef2SBarry Smith .  initial_time - the initial time
987d763cef2SBarry Smith -  time_step - the size of the timestep
988d763cef2SBarry Smith 
989d763cef2SBarry Smith    Level: intermediate
990d763cef2SBarry Smith 
991d763cef2SBarry Smith .seealso: TSSetTimeStep(), TSGetTimeStep()
992d763cef2SBarry Smith 
993d763cef2SBarry Smith .keywords: TS, set, initial, timestep
994d763cef2SBarry Smith @*/
9957087cfbeSBarry Smith PetscErrorCode  TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step)
996d763cef2SBarry Smith {
997e144a568SJed Brown   PetscErrorCode ierr;
998e144a568SJed Brown 
999d763cef2SBarry Smith   PetscFunctionBegin;
10000700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1001e144a568SJed Brown   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
1002d8cd7023SBarry Smith   ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr);
1003d763cef2SBarry Smith   PetscFunctionReturn(0);
1004d763cef2SBarry Smith }
1005d763cef2SBarry Smith 
10064a2ae208SSatish Balay #undef __FUNCT__
10074a2ae208SSatish Balay #define __FUNCT__ "TSSetTimeStep"
1008d763cef2SBarry Smith /*@
1009d763cef2SBarry Smith    TSSetTimeStep - Allows one to reset the timestep at any time,
1010d763cef2SBarry Smith    useful for simple pseudo-timestepping codes.
1011d763cef2SBarry Smith 
10123f9fe445SBarry Smith    Logically Collective on TS
1013d763cef2SBarry Smith 
1014d763cef2SBarry Smith    Input Parameters:
1015d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1016d763cef2SBarry Smith -  time_step - the size of the timestep
1017d763cef2SBarry Smith 
1018d763cef2SBarry Smith    Level: intermediate
1019d763cef2SBarry Smith 
1020d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1021d763cef2SBarry Smith 
1022d763cef2SBarry Smith .keywords: TS, set, timestep
1023d763cef2SBarry Smith @*/
10247087cfbeSBarry Smith PetscErrorCode  TSSetTimeStep(TS ts,PetscReal time_step)
1025d763cef2SBarry Smith {
1026d763cef2SBarry Smith   PetscFunctionBegin;
10270700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1028c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,time_step,2);
1029d763cef2SBarry Smith   ts->time_step = time_step;
1030d763cef2SBarry Smith   PetscFunctionReturn(0);
1031d763cef2SBarry Smith }
1032d763cef2SBarry Smith 
10334a2ae208SSatish Balay #undef __FUNCT__
1034a43b19c4SJed Brown #define __FUNCT__ "TSSetExactFinalTime"
1035a43b19c4SJed Brown /*@
1036a43b19c4SJed Brown    TSSetExactFinalTime - Determines whether to interpolate solution to the
1037a43b19c4SJed Brown       exact final time requested by the user or just returns it at the final time
1038a43b19c4SJed Brown       it computed.
1039a43b19c4SJed Brown 
1040a43b19c4SJed Brown   Logically Collective on TS
1041a43b19c4SJed Brown 
1042a43b19c4SJed Brown    Input Parameter:
1043a43b19c4SJed Brown +   ts - the time-step context
1044a43b19c4SJed Brown -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
1045a43b19c4SJed Brown 
1046a43b19c4SJed Brown    Level: beginner
1047a43b19c4SJed Brown 
1048a43b19c4SJed Brown .seealso: TSSetDuration()
1049a43b19c4SJed Brown @*/
1050a43b19c4SJed Brown PetscErrorCode  TSSetExactFinalTime(TS ts,PetscBool flg)
1051a43b19c4SJed Brown {
1052a43b19c4SJed Brown 
1053a43b19c4SJed Brown   PetscFunctionBegin;
1054a43b19c4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1055d8cd7023SBarry Smith   PetscValidLogicalCollectiveBool(ts,flg,2);
1056a43b19c4SJed Brown   ts->exact_final_time = flg;
1057a43b19c4SJed Brown   PetscFunctionReturn(0);
1058a43b19c4SJed Brown }
1059a43b19c4SJed Brown 
1060a43b19c4SJed Brown #undef __FUNCT__
10614a2ae208SSatish Balay #define __FUNCT__ "TSGetTimeStep"
1062d763cef2SBarry Smith /*@
1063d763cef2SBarry Smith    TSGetTimeStep - Gets the current timestep size.
1064d763cef2SBarry Smith 
1065d763cef2SBarry Smith    Not Collective
1066d763cef2SBarry Smith 
1067d763cef2SBarry Smith    Input Parameter:
1068d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1069d763cef2SBarry Smith 
1070d763cef2SBarry Smith    Output Parameter:
1071d763cef2SBarry Smith .  dt - the current timestep size
1072d763cef2SBarry Smith 
1073d763cef2SBarry Smith    Level: intermediate
1074d763cef2SBarry Smith 
1075d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
1076d763cef2SBarry Smith 
1077d763cef2SBarry Smith .keywords: TS, get, timestep
1078d763cef2SBarry Smith @*/
10797087cfbeSBarry Smith PetscErrorCode  TSGetTimeStep(TS ts,PetscReal* dt)
1080d763cef2SBarry Smith {
1081d763cef2SBarry Smith   PetscFunctionBegin;
10820700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10834482741eSBarry Smith   PetscValidDoublePointer(dt,2);
1084d763cef2SBarry Smith   *dt = ts->time_step;
1085d763cef2SBarry Smith   PetscFunctionReturn(0);
1086d763cef2SBarry Smith }
1087d763cef2SBarry Smith 
10884a2ae208SSatish Balay #undef __FUNCT__
10894a2ae208SSatish Balay #define __FUNCT__ "TSGetSolution"
1090d8e5e3e6SSatish Balay /*@
1091d763cef2SBarry Smith    TSGetSolution - Returns the solution at the present timestep. It
1092d763cef2SBarry Smith    is valid to call this routine inside the function that you are evaluating
1093d763cef2SBarry Smith    in order to move to the new timestep. This vector not changed until
1094d763cef2SBarry Smith    the solution at the next timestep has been calculated.
1095d763cef2SBarry Smith 
1096d763cef2SBarry Smith    Not Collective, but Vec returned is parallel if TS is parallel
1097d763cef2SBarry Smith 
1098d763cef2SBarry Smith    Input Parameter:
1099d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1100d763cef2SBarry Smith 
1101d763cef2SBarry Smith    Output Parameter:
1102d763cef2SBarry Smith .  v - the vector containing the solution
1103d763cef2SBarry Smith 
1104d763cef2SBarry Smith    Level: intermediate
1105d763cef2SBarry Smith 
1106d763cef2SBarry Smith .seealso: TSGetTimeStep()
1107d763cef2SBarry Smith 
1108d763cef2SBarry Smith .keywords: TS, timestep, get, solution
1109d763cef2SBarry Smith @*/
11107087cfbeSBarry Smith PetscErrorCode  TSGetSolution(TS ts,Vec *v)
1111d763cef2SBarry Smith {
1112d763cef2SBarry Smith   PetscFunctionBegin;
11130700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
11144482741eSBarry Smith   PetscValidPointer(v,2);
11158737fe31SLisandro Dalcin   *v = ts->vec_sol;
1116d763cef2SBarry Smith   PetscFunctionReturn(0);
1117d763cef2SBarry Smith }
1118d763cef2SBarry Smith 
1119bdad233fSMatthew Knepley /* ----- Routines to initialize and destroy a timestepper ---- */
11204a2ae208SSatish Balay #undef __FUNCT__
1121bdad233fSMatthew Knepley #define __FUNCT__ "TSSetProblemType"
1122d8e5e3e6SSatish Balay /*@
1123bdad233fSMatthew Knepley   TSSetProblemType - Sets the type of problem to be solved.
1124d763cef2SBarry Smith 
1125bdad233fSMatthew Knepley   Not collective
1126d763cef2SBarry Smith 
1127bdad233fSMatthew Knepley   Input Parameters:
1128bdad233fSMatthew Knepley + ts   - The TS
1129bdad233fSMatthew Knepley - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1130d763cef2SBarry Smith .vb
1131d763cef2SBarry Smith          U_t = A U
1132d763cef2SBarry Smith          U_t = A(t) U
1133d763cef2SBarry Smith          U_t = F(t,U)
1134d763cef2SBarry Smith .ve
1135d763cef2SBarry Smith 
1136d763cef2SBarry Smith    Level: beginner
1137d763cef2SBarry Smith 
1138bdad233fSMatthew Knepley .keywords: TS, problem type
1139bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1140d763cef2SBarry Smith @*/
11417087cfbeSBarry Smith PetscErrorCode  TSSetProblemType(TS ts, TSProblemType type)
1142a7cc72afSBarry Smith {
11439e2a6581SJed Brown   PetscErrorCode ierr;
11449e2a6581SJed Brown 
1145d763cef2SBarry Smith   PetscFunctionBegin;
11460700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1147bdad233fSMatthew Knepley   ts->problem_type = type;
11489e2a6581SJed Brown   if (type == TS_LINEAR) {
11499e2a6581SJed Brown     SNES snes;
11509e2a6581SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
11519e2a6581SJed Brown     ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
11529e2a6581SJed Brown   }
1153d763cef2SBarry Smith   PetscFunctionReturn(0);
1154d763cef2SBarry Smith }
1155d763cef2SBarry Smith 
1156bdad233fSMatthew Knepley #undef __FUNCT__
1157bdad233fSMatthew Knepley #define __FUNCT__ "TSGetProblemType"
1158bdad233fSMatthew Knepley /*@C
1159bdad233fSMatthew Knepley   TSGetProblemType - Gets the type of problem to be solved.
1160bdad233fSMatthew Knepley 
1161bdad233fSMatthew Knepley   Not collective
1162bdad233fSMatthew Knepley 
1163bdad233fSMatthew Knepley   Input Parameter:
1164bdad233fSMatthew Knepley . ts   - The TS
1165bdad233fSMatthew Knepley 
1166bdad233fSMatthew Knepley   Output Parameter:
1167bdad233fSMatthew Knepley . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms
1168bdad233fSMatthew Knepley .vb
1169089b2837SJed Brown          M U_t = A U
1170089b2837SJed Brown          M(t) U_t = A(t) U
1171bdad233fSMatthew Knepley          U_t = F(t,U)
1172bdad233fSMatthew Knepley .ve
1173bdad233fSMatthew Knepley 
1174bdad233fSMatthew Knepley    Level: beginner
1175bdad233fSMatthew Knepley 
1176bdad233fSMatthew Knepley .keywords: TS, problem type
1177bdad233fSMatthew Knepley .seealso: TSSetUp(), TSProblemType, TS
1178bdad233fSMatthew Knepley @*/
11797087cfbeSBarry Smith PetscErrorCode  TSGetProblemType(TS ts, TSProblemType *type)
1180a7cc72afSBarry Smith {
1181bdad233fSMatthew Knepley   PetscFunctionBegin;
11820700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
11834482741eSBarry Smith   PetscValidIntPointer(type,2);
1184bdad233fSMatthew Knepley   *type = ts->problem_type;
1185bdad233fSMatthew Knepley   PetscFunctionReturn(0);
1186bdad233fSMatthew Knepley }
1187d763cef2SBarry Smith 
11884a2ae208SSatish Balay #undef __FUNCT__
11894a2ae208SSatish Balay #define __FUNCT__ "TSSetUp"
1190d763cef2SBarry Smith /*@
1191d763cef2SBarry Smith    TSSetUp - Sets up the internal data structures for the later use
1192d763cef2SBarry Smith    of a timestepper.
1193d763cef2SBarry Smith 
1194d763cef2SBarry Smith    Collective on TS
1195d763cef2SBarry Smith 
1196d763cef2SBarry Smith    Input Parameter:
1197d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1198d763cef2SBarry Smith 
1199d763cef2SBarry Smith    Notes:
1200d763cef2SBarry Smith    For basic use of the TS solvers the user need not explicitly call
1201d763cef2SBarry Smith    TSSetUp(), since these actions will automatically occur during
1202d763cef2SBarry Smith    the call to TSStep().  However, if one wishes to control this
1203d763cef2SBarry Smith    phase separately, TSSetUp() should be called after TSCreate()
1204d763cef2SBarry Smith    and optional routines of the form TSSetXXX(), but before TSStep().
1205d763cef2SBarry Smith 
1206d763cef2SBarry Smith    Level: advanced
1207d763cef2SBarry Smith 
1208d763cef2SBarry Smith .keywords: TS, timestep, setup
1209d763cef2SBarry Smith 
1210d763cef2SBarry Smith .seealso: TSCreate(), TSStep(), TSDestroy()
1211d763cef2SBarry Smith @*/
12127087cfbeSBarry Smith PetscErrorCode  TSSetUp(TS ts)
1213d763cef2SBarry Smith {
1214dfbe8321SBarry Smith   PetscErrorCode ierr;
1215d763cef2SBarry Smith 
1216d763cef2SBarry Smith   PetscFunctionBegin;
12170700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1218277b19d0SLisandro Dalcin   if (ts->setupcalled) PetscFunctionReturn(0);
1219277b19d0SLisandro Dalcin 
12207adad957SLisandro Dalcin   if (!((PetscObject)ts)->type_name) {
12219596e0b4SJed Brown     ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
1222d763cef2SBarry Smith   }
1223a43b19c4SJed Brown   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1224277b19d0SLisandro Dalcin 
1225277b19d0SLisandro Dalcin   if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first");
1226277b19d0SLisandro Dalcin 
12271c3436cfSJed Brown   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
12281c3436cfSJed Brown 
1229277b19d0SLisandro Dalcin   if (ts->ops->setup) {
1230000e7ae3SMatthew Knepley     ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr);
1231277b19d0SLisandro Dalcin   }
1232277b19d0SLisandro Dalcin 
1233277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_TRUE;
1234277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1235277b19d0SLisandro Dalcin }
1236277b19d0SLisandro Dalcin 
1237277b19d0SLisandro Dalcin #undef __FUNCT__
1238277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset"
1239277b19d0SLisandro Dalcin /*@
1240277b19d0SLisandro Dalcin    TSReset - Resets a TS context and removes any allocated Vecs and Mats.
1241277b19d0SLisandro Dalcin 
1242277b19d0SLisandro Dalcin    Collective on TS
1243277b19d0SLisandro Dalcin 
1244277b19d0SLisandro Dalcin    Input Parameter:
1245277b19d0SLisandro Dalcin .  ts - the TS context obtained from TSCreate()
1246277b19d0SLisandro Dalcin 
1247277b19d0SLisandro Dalcin    Level: beginner
1248277b19d0SLisandro Dalcin 
1249277b19d0SLisandro Dalcin .keywords: TS, timestep, reset
1250277b19d0SLisandro Dalcin 
1251277b19d0SLisandro Dalcin .seealso: TSCreate(), TSSetup(), TSDestroy()
1252277b19d0SLisandro Dalcin @*/
1253277b19d0SLisandro Dalcin PetscErrorCode  TSReset(TS ts)
1254277b19d0SLisandro Dalcin {
1255277b19d0SLisandro Dalcin   PetscErrorCode ierr;
1256277b19d0SLisandro Dalcin 
1257277b19d0SLisandro Dalcin   PetscFunctionBegin;
1258277b19d0SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1259277b19d0SLisandro Dalcin   if (ts->ops->reset) {
1260277b19d0SLisandro Dalcin     ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr);
1261277b19d0SLisandro Dalcin   }
1262277b19d0SLisandro Dalcin   if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);}
12634e684422SJed Brown   ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr);
12644e684422SJed Brown   ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr);
1265214bc6a2SJed Brown   ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr);
12666bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
1267e3d84a46SJed Brown   ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
1268e3d84a46SJed Brown   ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
126938637c2eSJed Brown   ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);
1270277b19d0SLisandro Dalcin   ts->setupcalled = PETSC_FALSE;
1271d763cef2SBarry Smith   PetscFunctionReturn(0);
1272d763cef2SBarry Smith }
1273d763cef2SBarry Smith 
12744a2ae208SSatish Balay #undef __FUNCT__
12754a2ae208SSatish Balay #define __FUNCT__ "TSDestroy"
1276d8e5e3e6SSatish Balay /*@
1277d763cef2SBarry Smith    TSDestroy - Destroys the timestepper context that was created
1278d763cef2SBarry Smith    with TSCreate().
1279d763cef2SBarry Smith 
1280d763cef2SBarry Smith    Collective on TS
1281d763cef2SBarry Smith 
1282d763cef2SBarry Smith    Input Parameter:
1283d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1284d763cef2SBarry Smith 
1285d763cef2SBarry Smith    Level: beginner
1286d763cef2SBarry Smith 
1287d763cef2SBarry Smith .keywords: TS, timestepper, destroy
1288d763cef2SBarry Smith 
1289d763cef2SBarry Smith .seealso: TSCreate(), TSSetUp(), TSSolve()
1290d763cef2SBarry Smith @*/
12916bf464f9SBarry Smith PetscErrorCode  TSDestroy(TS *ts)
1292d763cef2SBarry Smith {
12936849ba73SBarry Smith   PetscErrorCode ierr;
1294d763cef2SBarry Smith 
1295d763cef2SBarry Smith   PetscFunctionBegin;
12966bf464f9SBarry Smith   if (!*ts) PetscFunctionReturn(0);
12976bf464f9SBarry Smith   PetscValidHeaderSpecific((*ts),TS_CLASSID,1);
12986bf464f9SBarry Smith   if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);}
1299d763cef2SBarry Smith 
13006bf464f9SBarry Smith   ierr = TSReset((*ts));CHKERRQ(ierr);
1301277b19d0SLisandro Dalcin 
1302be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
13036bf464f9SBarry Smith   ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr);
13046bf464f9SBarry Smith   if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);}
13056d4c513bSLisandro Dalcin 
130684df9cb4SJed Brown   ierr = TSAdaptDestroy(&(*ts)->adapt);CHKERRQ(ierr);
13076bf464f9SBarry Smith   ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr);
13086bf464f9SBarry Smith   ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr);
13096bf464f9SBarry Smith   ierr = TSMonitorCancel((*ts));CHKERRQ(ierr);
13106d4c513bSLisandro Dalcin 
13113b5f76d0SSean Farley   ierr = PetscFree((*ts)->userops);
13123b5f76d0SSean Farley 
1313a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr);
1314d763cef2SBarry Smith   PetscFunctionReturn(0);
1315d763cef2SBarry Smith }
1316d763cef2SBarry Smith 
13174a2ae208SSatish Balay #undef __FUNCT__
13184a2ae208SSatish Balay #define __FUNCT__ "TSGetSNES"
1319d8e5e3e6SSatish Balay /*@
1320d763cef2SBarry Smith    TSGetSNES - Returns the SNES (nonlinear solver) associated with
1321d763cef2SBarry Smith    a TS (timestepper) context. Valid only for nonlinear problems.
1322d763cef2SBarry Smith 
1323d763cef2SBarry Smith    Not Collective, but SNES is parallel if TS is parallel
1324d763cef2SBarry Smith 
1325d763cef2SBarry Smith    Input Parameter:
1326d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1327d763cef2SBarry Smith 
1328d763cef2SBarry Smith    Output Parameter:
1329d763cef2SBarry Smith .  snes - the nonlinear solver context
1330d763cef2SBarry Smith 
1331d763cef2SBarry Smith    Notes:
1332d763cef2SBarry Smith    The user can then directly manipulate the SNES context to set various
1333d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
133494b7f48cSBarry Smith    KSP, KSP, and PC contexts as well.
1335d763cef2SBarry Smith 
1336d763cef2SBarry Smith    TSGetSNES() does not work for integrators that do not use SNES; in
1337d763cef2SBarry Smith    this case TSGetSNES() returns PETSC_NULL in snes.
1338d763cef2SBarry Smith 
1339d763cef2SBarry Smith    Level: beginner
1340d763cef2SBarry Smith 
1341d763cef2SBarry Smith .keywords: timestep, get, SNES
1342d763cef2SBarry Smith @*/
13437087cfbeSBarry Smith PetscErrorCode  TSGetSNES(TS ts,SNES *snes)
1344d763cef2SBarry Smith {
1345d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1346d372ba47SLisandro Dalcin 
1347d763cef2SBarry Smith   PetscFunctionBegin;
13480700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13494482741eSBarry Smith   PetscValidPointer(snes,2);
1350d372ba47SLisandro Dalcin   if (!ts->snes) {
1351d372ba47SLisandro Dalcin     ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr);
1352d372ba47SLisandro Dalcin     ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr);
1353d372ba47SLisandro Dalcin     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr);
1354496e6a7aSJed Brown     if (ts->dm) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
13559e2a6581SJed Brown     if (ts->problem_type == TS_LINEAR) {
13569e2a6581SJed Brown       ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr);
13579e2a6581SJed Brown     }
1358d372ba47SLisandro Dalcin   }
1359d763cef2SBarry Smith   *snes = ts->snes;
1360d763cef2SBarry Smith   PetscFunctionReturn(0);
1361d763cef2SBarry Smith }
1362d763cef2SBarry Smith 
13634a2ae208SSatish Balay #undef __FUNCT__
136494b7f48cSBarry Smith #define __FUNCT__ "TSGetKSP"
1365d8e5e3e6SSatish Balay /*@
136694b7f48cSBarry Smith    TSGetKSP - Returns the KSP (linear solver) associated with
1367d763cef2SBarry Smith    a TS (timestepper) context.
1368d763cef2SBarry Smith 
136994b7f48cSBarry Smith    Not Collective, but KSP is parallel if TS is parallel
1370d763cef2SBarry Smith 
1371d763cef2SBarry Smith    Input Parameter:
1372d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1373d763cef2SBarry Smith 
1374d763cef2SBarry Smith    Output Parameter:
137594b7f48cSBarry Smith .  ksp - the nonlinear solver context
1376d763cef2SBarry Smith 
1377d763cef2SBarry Smith    Notes:
137894b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
1379d763cef2SBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
1380d763cef2SBarry Smith    KSP and PC contexts as well.
1381d763cef2SBarry Smith 
138294b7f48cSBarry Smith    TSGetKSP() does not work for integrators that do not use KSP;
138394b7f48cSBarry Smith    in this case TSGetKSP() returns PETSC_NULL in ksp.
1384d763cef2SBarry Smith 
1385d763cef2SBarry Smith    Level: beginner
1386d763cef2SBarry Smith 
138794b7f48cSBarry Smith .keywords: timestep, get, KSP
1388d763cef2SBarry Smith @*/
13897087cfbeSBarry Smith PetscErrorCode  TSGetKSP(TS ts,KSP *ksp)
1390d763cef2SBarry Smith {
1391d372ba47SLisandro Dalcin   PetscErrorCode ierr;
1392089b2837SJed Brown   SNES           snes;
1393d372ba47SLisandro Dalcin 
1394d763cef2SBarry Smith   PetscFunctionBegin;
13950700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13964482741eSBarry Smith   PetscValidPointer(ksp,2);
139717186662SBarry Smith   if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first");
1398e32f2f54SBarry Smith   if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()");
1399089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1400089b2837SJed Brown   ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr);
1401d763cef2SBarry Smith   PetscFunctionReturn(0);
1402d763cef2SBarry Smith }
1403d763cef2SBarry Smith 
1404d763cef2SBarry Smith /* ----------- Routines to set solver parameters ---------- */
1405d763cef2SBarry Smith 
14064a2ae208SSatish Balay #undef __FUNCT__
1407adb62b0dSMatthew Knepley #define __FUNCT__ "TSGetDuration"
1408adb62b0dSMatthew Knepley /*@
1409adb62b0dSMatthew Knepley    TSGetDuration - Gets the maximum number of timesteps to use and
1410adb62b0dSMatthew Knepley    maximum time for iteration.
1411adb62b0dSMatthew Knepley 
14123f9fe445SBarry Smith    Not Collective
1413adb62b0dSMatthew Knepley 
1414adb62b0dSMatthew Knepley    Input Parameters:
1415adb62b0dSMatthew Knepley +  ts       - the TS context obtained from TSCreate()
1416adb62b0dSMatthew Knepley .  maxsteps - maximum number of iterations to use, or PETSC_NULL
1417adb62b0dSMatthew Knepley -  maxtime  - final time to iterate to, or PETSC_NULL
1418adb62b0dSMatthew Knepley 
1419adb62b0dSMatthew Knepley    Level: intermediate
1420adb62b0dSMatthew Knepley 
1421adb62b0dSMatthew Knepley .keywords: TS, timestep, get, maximum, iterations, time
1422adb62b0dSMatthew Knepley @*/
14237087cfbeSBarry Smith PetscErrorCode  TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime)
1424adb62b0dSMatthew Knepley {
1425adb62b0dSMatthew Knepley   PetscFunctionBegin;
14260700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1427abc0a331SBarry Smith   if (maxsteps) {
14284482741eSBarry Smith     PetscValidIntPointer(maxsteps,2);
1429adb62b0dSMatthew Knepley     *maxsteps = ts->max_steps;
1430adb62b0dSMatthew Knepley   }
1431abc0a331SBarry Smith   if (maxtime) {
14324482741eSBarry Smith     PetscValidScalarPointer(maxtime,3);
1433adb62b0dSMatthew Knepley     *maxtime  = ts->max_time;
1434adb62b0dSMatthew Knepley   }
1435adb62b0dSMatthew Knepley   PetscFunctionReturn(0);
1436adb62b0dSMatthew Knepley }
1437adb62b0dSMatthew Knepley 
1438adb62b0dSMatthew Knepley #undef __FUNCT__
14394a2ae208SSatish Balay #define __FUNCT__ "TSSetDuration"
1440d763cef2SBarry Smith /*@
1441d763cef2SBarry Smith    TSSetDuration - Sets the maximum number of timesteps to use and
1442d763cef2SBarry Smith    maximum time for iteration.
1443d763cef2SBarry Smith 
14443f9fe445SBarry Smith    Logically Collective on TS
1445d763cef2SBarry Smith 
1446d763cef2SBarry Smith    Input Parameters:
1447d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1448d763cef2SBarry Smith .  maxsteps - maximum number of iterations to use
1449d763cef2SBarry Smith -  maxtime - final time to iterate to
1450d763cef2SBarry Smith 
1451d763cef2SBarry Smith    Options Database Keys:
1452d763cef2SBarry Smith .  -ts_max_steps <maxsteps> - Sets maxsteps
14533bca7d26SBarry Smith .  -ts_final_time <maxtime> - Sets maxtime
1454d763cef2SBarry Smith 
1455d763cef2SBarry Smith    Notes:
1456d763cef2SBarry Smith    The default maximum number of iterations is 5000. Default time is 5.0
1457d763cef2SBarry Smith 
1458d763cef2SBarry Smith    Level: intermediate
1459d763cef2SBarry Smith 
1460d763cef2SBarry Smith .keywords: TS, timestep, set, maximum, iterations
1461a43b19c4SJed Brown 
1462a43b19c4SJed Brown .seealso: TSSetExactFinalTime()
1463d763cef2SBarry Smith @*/
14647087cfbeSBarry Smith PetscErrorCode  TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime)
1465d763cef2SBarry Smith {
1466d763cef2SBarry Smith   PetscFunctionBegin;
14670700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1468c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,maxsteps,2);
1469c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,maxtime,2);
147039b7ec4bSSean Farley   if (maxsteps >= 0) ts->max_steps = maxsteps;
147139b7ec4bSSean Farley   if (maxtime != PETSC_DEFAULT) ts->max_time  = maxtime;
1472d763cef2SBarry Smith   PetscFunctionReturn(0);
1473d763cef2SBarry Smith }
1474d763cef2SBarry Smith 
14754a2ae208SSatish Balay #undef __FUNCT__
14764a2ae208SSatish Balay #define __FUNCT__ "TSSetSolution"
1477d763cef2SBarry Smith /*@
1478d763cef2SBarry Smith    TSSetSolution - Sets the initial solution vector
1479d763cef2SBarry Smith    for use by the TS routines.
1480d763cef2SBarry Smith 
14813f9fe445SBarry Smith    Logically Collective on TS and Vec
1482d763cef2SBarry Smith 
1483d763cef2SBarry Smith    Input Parameters:
1484d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1485d763cef2SBarry Smith -  x - the solution vector
1486d763cef2SBarry Smith 
1487d763cef2SBarry Smith    Level: beginner
1488d763cef2SBarry Smith 
1489d763cef2SBarry Smith .keywords: TS, timestep, set, solution, initial conditions
1490d763cef2SBarry Smith @*/
14917087cfbeSBarry Smith PetscErrorCode  TSSetSolution(TS ts,Vec x)
1492d763cef2SBarry Smith {
14938737fe31SLisandro Dalcin   PetscErrorCode ierr;
1494496e6a7aSJed Brown   DM             dm;
14958737fe31SLisandro Dalcin 
1496d763cef2SBarry Smith   PetscFunctionBegin;
14970700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
14980700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
14998737fe31SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
15006bf464f9SBarry Smith   ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr);
15018737fe31SLisandro Dalcin   ts->vec_sol = x;
1502496e6a7aSJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1503496e6a7aSJed Brown   ierr = DMShellSetGlobalVector(dm,x);CHKERRQ(ierr);
1504d763cef2SBarry Smith   PetscFunctionReturn(0);
1505d763cef2SBarry Smith }
1506d763cef2SBarry Smith 
1507e74ef692SMatthew Knepley #undef __FUNCT__
1508e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPreStep"
1509ac226902SBarry Smith /*@C
1510000e7ae3SMatthew Knepley   TSSetPreStep - Sets the general-purpose function
15113f2090d5SJed Brown   called once at the beginning of each time step.
1512000e7ae3SMatthew Knepley 
15133f9fe445SBarry Smith   Logically Collective on TS
1514000e7ae3SMatthew Knepley 
1515000e7ae3SMatthew Knepley   Input Parameters:
1516000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1517000e7ae3SMatthew Knepley - func - The function
1518000e7ae3SMatthew Knepley 
1519000e7ae3SMatthew Knepley   Calling sequence of func:
1520000e7ae3SMatthew Knepley . func (TS ts);
1521000e7ae3SMatthew Knepley 
1522000e7ae3SMatthew Knepley   Level: intermediate
1523000e7ae3SMatthew Knepley 
1524*b8123daeSJed Brown   Note:
1525*b8123daeSJed Brown   If a step is rejected, TSStep() will call this routine again before each attempt.
1526*b8123daeSJed Brown   The last completed time step number can be queried using TSGetTimeStepNumber(), the
1527*b8123daeSJed Brown   size of the step being attempted can be obtained using TSGetTimeStep().
1528*b8123daeSJed Brown 
1529000e7ae3SMatthew Knepley .keywords: TS, timestep
1530*b8123daeSJed Brown .seealso: TSSetPreStage(), TSSetPostStep(), TSStep()
1531000e7ae3SMatthew Knepley @*/
15327087cfbeSBarry Smith PetscErrorCode  TSSetPreStep(TS ts, PetscErrorCode (*func)(TS))
1533000e7ae3SMatthew Knepley {
1534000e7ae3SMatthew Knepley   PetscFunctionBegin;
15350700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1536000e7ae3SMatthew Knepley   ts->ops->prestep = func;
1537000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1538000e7ae3SMatthew Knepley }
1539000e7ae3SMatthew Knepley 
1540e74ef692SMatthew Knepley #undef __FUNCT__
15413f2090d5SJed Brown #define __FUNCT__ "TSPreStep"
154209ee8438SJed Brown /*@
15433f2090d5SJed Brown   TSPreStep - Runs the user-defined pre-step function.
15443f2090d5SJed Brown 
15453f2090d5SJed Brown   Collective on TS
15463f2090d5SJed Brown 
15473f2090d5SJed Brown   Input Parameters:
15483f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
15493f2090d5SJed Brown 
15503f2090d5SJed Brown   Notes:
15513f2090d5SJed Brown   TSPreStep() is typically used within time stepping implementations,
15523f2090d5SJed Brown   so most users would not generally call this routine themselves.
15533f2090d5SJed Brown 
15543f2090d5SJed Brown   Level: developer
15553f2090d5SJed Brown 
15563f2090d5SJed Brown .keywords: TS, timestep
1557*b8123daeSJed Brown .seealso: TSSetPreStep(), TSPreStage(), TSPostStep()
15583f2090d5SJed Brown @*/
15597087cfbeSBarry Smith PetscErrorCode  TSPreStep(TS ts)
15603f2090d5SJed Brown {
15613f2090d5SJed Brown   PetscErrorCode ierr;
15623f2090d5SJed Brown 
15633f2090d5SJed Brown   PetscFunctionBegin;
15640700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
156572ac3e02SJed Brown   if (ts->ops->prestep) {
15663f2090d5SJed Brown     PetscStackPush("TS PreStep function");
15673f2090d5SJed Brown     ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr);
15683f2090d5SJed Brown     PetscStackPop;
1569312ce896SJed Brown   }
15703f2090d5SJed Brown   PetscFunctionReturn(0);
15713f2090d5SJed Brown }
15723f2090d5SJed Brown 
15733f2090d5SJed Brown #undef __FUNCT__
1574*b8123daeSJed Brown #define __FUNCT__ "TSSetPreStage"
1575*b8123daeSJed Brown /*@C
1576*b8123daeSJed Brown   TSSetPreStage - Sets the general-purpose function
1577*b8123daeSJed Brown   called once at the beginning of each stage.
1578*b8123daeSJed Brown 
1579*b8123daeSJed Brown   Logically Collective on TS
1580*b8123daeSJed Brown 
1581*b8123daeSJed Brown   Input Parameters:
1582*b8123daeSJed Brown + ts   - The TS context obtained from TSCreate()
1583*b8123daeSJed Brown - func - The function
1584*b8123daeSJed Brown 
1585*b8123daeSJed Brown   Calling sequence of func:
1586*b8123daeSJed Brown . PetscErrorCode func(TS ts, PetscReal stagetime);
1587*b8123daeSJed Brown 
1588*b8123daeSJed Brown   Level: intermediate
1589*b8123daeSJed Brown 
1590*b8123daeSJed Brown   Note:
1591*b8123daeSJed Brown   There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried.
1592*b8123daeSJed Brown   The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being
1593*b8123daeSJed Brown   attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime().
1594*b8123daeSJed Brown 
1595*b8123daeSJed Brown .keywords: TS, timestep
1596*b8123daeSJed Brown .seealso: TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext()
1597*b8123daeSJed Brown @*/
1598*b8123daeSJed Brown PetscErrorCode  TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal))
1599*b8123daeSJed Brown {
1600*b8123daeSJed Brown   PetscFunctionBegin;
1601*b8123daeSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1602*b8123daeSJed Brown   ts->ops->prestage = func;
1603*b8123daeSJed Brown   PetscFunctionReturn(0);
1604*b8123daeSJed Brown }
1605*b8123daeSJed Brown 
1606*b8123daeSJed Brown #undef __FUNCT__
1607*b8123daeSJed Brown #define __FUNCT__ "TSPreStage"
1608*b8123daeSJed Brown /*@
1609*b8123daeSJed Brown   TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage()
1610*b8123daeSJed Brown 
1611*b8123daeSJed Brown   Collective on TS
1612*b8123daeSJed Brown 
1613*b8123daeSJed Brown   Input Parameters:
1614*b8123daeSJed Brown . ts   - The TS context obtained from TSCreate()
1615*b8123daeSJed Brown 
1616*b8123daeSJed Brown   Notes:
1617*b8123daeSJed Brown   TSPreStage() is typically used within time stepping implementations,
1618*b8123daeSJed Brown   most users would not generally call this routine themselves.
1619*b8123daeSJed Brown 
1620*b8123daeSJed Brown   Level: developer
1621*b8123daeSJed Brown 
1622*b8123daeSJed Brown .keywords: TS, timestep
1623*b8123daeSJed Brown .seealso: TSSetPreStep(), TSPreStep(), TSPostStep()
1624*b8123daeSJed Brown @*/
1625*b8123daeSJed Brown PetscErrorCode  TSPreStage(TS ts, PetscReal stagetime)
1626*b8123daeSJed Brown {
1627*b8123daeSJed Brown   PetscErrorCode ierr;
1628*b8123daeSJed Brown 
1629*b8123daeSJed Brown   PetscFunctionBegin;
1630*b8123daeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1631*b8123daeSJed Brown   if (ts->ops->prestage) {
1632*b8123daeSJed Brown     PetscStackPush("TS PreStage function");
1633*b8123daeSJed Brown     ierr = (*ts->ops->prestage)(ts,stagetime);CHKERRQ(ierr);
1634*b8123daeSJed Brown     PetscStackPop;
1635*b8123daeSJed Brown   }
1636*b8123daeSJed Brown   PetscFunctionReturn(0);
1637*b8123daeSJed Brown }
1638*b8123daeSJed Brown 
1639*b8123daeSJed Brown #undef __FUNCT__
1640e74ef692SMatthew Knepley #define __FUNCT__ "TSSetPostStep"
1641ac226902SBarry Smith /*@C
1642000e7ae3SMatthew Knepley   TSSetPostStep - Sets the general-purpose function
16433f2090d5SJed Brown   called once at the end of each time step.
1644000e7ae3SMatthew Knepley 
16453f9fe445SBarry Smith   Logically Collective on TS
1646000e7ae3SMatthew Knepley 
1647000e7ae3SMatthew Knepley   Input Parameters:
1648000e7ae3SMatthew Knepley + ts   - The TS context obtained from TSCreate()
1649000e7ae3SMatthew Knepley - func - The function
1650000e7ae3SMatthew Knepley 
1651000e7ae3SMatthew Knepley   Calling sequence of func:
1652*b8123daeSJed Brown $ func (TS ts);
1653000e7ae3SMatthew Knepley 
1654000e7ae3SMatthew Knepley   Level: intermediate
1655000e7ae3SMatthew Knepley 
1656000e7ae3SMatthew Knepley .keywords: TS, timestep
1657*b8123daeSJed Brown .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime()
1658000e7ae3SMatthew Knepley @*/
16597087cfbeSBarry Smith PetscErrorCode  TSSetPostStep(TS ts, PetscErrorCode (*func)(TS))
1660000e7ae3SMatthew Knepley {
1661000e7ae3SMatthew Knepley   PetscFunctionBegin;
16620700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1663000e7ae3SMatthew Knepley   ts->ops->poststep = func;
1664000e7ae3SMatthew Knepley   PetscFunctionReturn(0);
1665000e7ae3SMatthew Knepley }
1666000e7ae3SMatthew Knepley 
1667e74ef692SMatthew Knepley #undef __FUNCT__
16683f2090d5SJed Brown #define __FUNCT__ "TSPostStep"
166909ee8438SJed Brown /*@
16703f2090d5SJed Brown   TSPostStep - Runs the user-defined post-step function.
16713f2090d5SJed Brown 
16723f2090d5SJed Brown   Collective on TS
16733f2090d5SJed Brown 
16743f2090d5SJed Brown   Input Parameters:
16753f2090d5SJed Brown . ts   - The TS context obtained from TSCreate()
16763f2090d5SJed Brown 
16773f2090d5SJed Brown   Notes:
16783f2090d5SJed Brown   TSPostStep() is typically used within time stepping implementations,
16793f2090d5SJed Brown   so most users would not generally call this routine themselves.
16803f2090d5SJed Brown 
16813f2090d5SJed Brown   Level: developer
16823f2090d5SJed Brown 
16833f2090d5SJed Brown .keywords: TS, timestep
16843f2090d5SJed Brown @*/
16857087cfbeSBarry Smith PetscErrorCode  TSPostStep(TS ts)
16863f2090d5SJed Brown {
16873f2090d5SJed Brown   PetscErrorCode ierr;
16883f2090d5SJed Brown 
16893f2090d5SJed Brown   PetscFunctionBegin;
16900700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
169172ac3e02SJed Brown   if (ts->ops->poststep) {
16923f2090d5SJed Brown     PetscStackPush("TS PostStep function");
16933f2090d5SJed Brown     ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr);
16943f2090d5SJed Brown     PetscStackPop;
169572ac3e02SJed Brown   }
16963f2090d5SJed Brown   PetscFunctionReturn(0);
16973f2090d5SJed Brown }
16983f2090d5SJed Brown 
1699d763cef2SBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
1700d763cef2SBarry Smith 
17014a2ae208SSatish Balay #undef __FUNCT__
1702a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSet"
1703d763cef2SBarry Smith /*@C
1704a6570f20SBarry Smith    TSMonitorSet - Sets an ADDITIONAL function that is to be used at every
1705d763cef2SBarry Smith    timestep to display the iteration's  progress.
1706d763cef2SBarry Smith 
17073f9fe445SBarry Smith    Logically Collective on TS
1708d763cef2SBarry Smith 
1709d763cef2SBarry Smith    Input Parameters:
1710d763cef2SBarry Smith +  ts - the TS context obtained from TSCreate()
1711e213d8f1SJed Brown .  monitor - monitoring routine
1712329f5518SBarry Smith .  mctx - [optional] user-defined context for private data for the
1713b3006f0bSLois Curfman McInnes              monitor routine (use PETSC_NULL if no context is desired)
1714b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1715b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
1716d763cef2SBarry Smith 
1717e213d8f1SJed Brown    Calling sequence of monitor:
1718e213d8f1SJed Brown $    int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx)
1719d763cef2SBarry Smith 
1720d763cef2SBarry Smith +    ts - the TS context
1721d763cef2SBarry Smith .    steps - iteration number
17221f06c33eSBarry Smith .    time - current time
1723d763cef2SBarry Smith .    x - current iterate
1724d763cef2SBarry Smith -    mctx - [optional] monitoring context
1725d763cef2SBarry Smith 
1726d763cef2SBarry Smith    Notes:
1727d763cef2SBarry Smith    This routine adds an additional monitor to the list of monitors that
1728d763cef2SBarry Smith    already has been loaded.
1729d763cef2SBarry Smith 
1730025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each TS object
1731025f1a04SBarry Smith 
1732d763cef2SBarry Smith    Level: intermediate
1733d763cef2SBarry Smith 
1734d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1735d763cef2SBarry Smith 
1736a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorCancel()
1737d763cef2SBarry Smith @*/
1738c2efdce3SBarry Smith PetscErrorCode  TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**))
1739d763cef2SBarry Smith {
1740d763cef2SBarry Smith   PetscFunctionBegin;
17410700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
174217186662SBarry Smith   if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1743d763cef2SBarry Smith   ts->monitor[ts->numbermonitors]           = monitor;
1744329f5518SBarry Smith   ts->mdestroy[ts->numbermonitors]          = mdestroy;
1745d763cef2SBarry Smith   ts->monitorcontext[ts->numbermonitors++]  = (void*)mctx;
1746d763cef2SBarry Smith   PetscFunctionReturn(0);
1747d763cef2SBarry Smith }
1748d763cef2SBarry Smith 
17494a2ae208SSatish Balay #undef __FUNCT__
1750a6570f20SBarry Smith #define __FUNCT__ "TSMonitorCancel"
1751d763cef2SBarry Smith /*@C
1752a6570f20SBarry Smith    TSMonitorCancel - Clears all the monitors that have been set on a time-step object.
1753d763cef2SBarry Smith 
17543f9fe445SBarry Smith    Logically Collective on TS
1755d763cef2SBarry Smith 
1756d763cef2SBarry Smith    Input Parameters:
1757d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1758d763cef2SBarry Smith 
1759d763cef2SBarry Smith    Notes:
1760d763cef2SBarry Smith    There is no way to remove a single, specific monitor.
1761d763cef2SBarry Smith 
1762d763cef2SBarry Smith    Level: intermediate
1763d763cef2SBarry Smith 
1764d763cef2SBarry Smith .keywords: TS, timestep, set, monitor
1765d763cef2SBarry Smith 
1766a6570f20SBarry Smith .seealso: TSMonitorDefault(), TSMonitorSet()
1767d763cef2SBarry Smith @*/
17687087cfbeSBarry Smith PetscErrorCode  TSMonitorCancel(TS ts)
1769d763cef2SBarry Smith {
1770d952e501SBarry Smith   PetscErrorCode ierr;
1771d952e501SBarry Smith   PetscInt       i;
1772d952e501SBarry Smith 
1773d763cef2SBarry Smith   PetscFunctionBegin;
17740700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1775d952e501SBarry Smith   for (i=0; i<ts->numbermonitors; i++) {
1776d952e501SBarry Smith     if (ts->mdestroy[i]) {
17773c4aec1bSBarry Smith       ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr);
1778d952e501SBarry Smith     }
1779d952e501SBarry Smith   }
1780d763cef2SBarry Smith   ts->numbermonitors = 0;
1781d763cef2SBarry Smith   PetscFunctionReturn(0);
1782d763cef2SBarry Smith }
1783d763cef2SBarry Smith 
17844a2ae208SSatish Balay #undef __FUNCT__
1785a6570f20SBarry Smith #define __FUNCT__ "TSMonitorDefault"
1786d8e5e3e6SSatish Balay /*@
1787a6570f20SBarry Smith    TSMonitorDefault - Sets the Default monitor
17885516499fSSatish Balay 
17895516499fSSatish Balay    Level: intermediate
179041251cbbSSatish Balay 
17915516499fSSatish Balay .keywords: TS, set, monitor
17925516499fSSatish Balay 
179341251cbbSSatish Balay .seealso: TSMonitorDefault(), TSMonitorSet()
179441251cbbSSatish Balay @*/
1795649052a6SBarry Smith PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy)
1796d763cef2SBarry Smith {
1797dfbe8321SBarry Smith   PetscErrorCode ierr;
1798649052a6SBarry Smith   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm);
1799d132466eSBarry Smith 
1800d763cef2SBarry Smith   PetscFunctionBegin;
1801649052a6SBarry Smith   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1802971832b6SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr);
1803649052a6SBarry Smith   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1804d763cef2SBarry Smith   PetscFunctionReturn(0);
1805d763cef2SBarry Smith }
1806d763cef2SBarry Smith 
18074a2ae208SSatish Balay #undef __FUNCT__
1808cd652676SJed Brown #define __FUNCT__ "TSSetRetainStages"
1809cd652676SJed Brown /*@
1810cd652676SJed Brown    TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available.
1811cd652676SJed Brown 
1812cd652676SJed Brown    Logically Collective on TS
1813cd652676SJed Brown 
1814cd652676SJed Brown    Input Argument:
1815cd652676SJed Brown .  ts - time stepping context
1816cd652676SJed Brown 
1817cd652676SJed Brown    Output Argument:
1818cd652676SJed Brown .  flg - PETSC_TRUE or PETSC_FALSE
1819cd652676SJed Brown 
1820cd652676SJed Brown    Level: intermediate
1821cd652676SJed Brown 
1822cd652676SJed Brown .keywords: TS, set
1823cd652676SJed Brown 
1824cd652676SJed Brown .seealso: TSInterpolate(), TSSetPostStep()
1825cd652676SJed Brown @*/
1826cd652676SJed Brown PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg)
1827cd652676SJed Brown {
1828cd652676SJed Brown 
1829cd652676SJed Brown   PetscFunctionBegin;
1830cd652676SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1831cd652676SJed Brown   ts->retain_stages = flg;
1832cd652676SJed Brown   PetscFunctionReturn(0);
1833cd652676SJed Brown }
1834cd652676SJed Brown 
1835cd652676SJed Brown #undef __FUNCT__
1836cd652676SJed Brown #define __FUNCT__ "TSInterpolate"
1837cd652676SJed Brown /*@
1838cd652676SJed Brown    TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval
1839cd652676SJed Brown 
1840cd652676SJed Brown    Collective on TS
1841cd652676SJed Brown 
1842cd652676SJed Brown    Input Argument:
1843cd652676SJed Brown +  ts - time stepping context
1844cd652676SJed Brown -  t - time to interpolate to
1845cd652676SJed Brown 
1846cd652676SJed Brown    Output Argument:
1847cd652676SJed Brown .  X - state at given time
1848cd652676SJed Brown 
1849cd652676SJed Brown    Notes:
1850cd652676SJed Brown    The user should call TSSetRetainStages() before taking a step in which interpolation will be requested.
1851cd652676SJed Brown 
1852cd652676SJed Brown    Level: intermediate
1853cd652676SJed Brown 
1854cd652676SJed Brown    Developer Notes:
1855cd652676SJed Brown    TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints.
1856cd652676SJed Brown 
1857cd652676SJed Brown .keywords: TS, set
1858cd652676SJed Brown 
1859cd652676SJed Brown .seealso: TSSetRetainStages(), TSSetPostStep()
1860cd652676SJed Brown @*/
1861cd652676SJed Brown PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X)
1862cd652676SJed Brown {
1863cd652676SJed Brown   PetscErrorCode ierr;
1864cd652676SJed Brown 
1865cd652676SJed Brown   PetscFunctionBegin;
1866cd652676SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1867d8cd7023SBarry 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);
1868cd652676SJed Brown   if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name);
1869cd652676SJed Brown   ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr);
1870cd652676SJed Brown   PetscFunctionReturn(0);
1871cd652676SJed Brown }
1872cd652676SJed Brown 
1873cd652676SJed Brown #undef __FUNCT__
18744a2ae208SSatish Balay #define __FUNCT__ "TSStep"
1875d763cef2SBarry Smith /*@
18766d9e5789SSean Farley    TSStep - Steps one time step
1877d763cef2SBarry Smith 
1878d763cef2SBarry Smith    Collective on TS
1879d763cef2SBarry Smith 
1880d763cef2SBarry Smith    Input Parameter:
1881d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
1882d763cef2SBarry Smith 
18836d9e5789SSean Farley    Level: intermediate
1884d763cef2SBarry Smith 
1885*b8123daeSJed Brown    Notes:
1886*b8123daeSJed Brown    The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may
1887*b8123daeSJed Brown    be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages.
1888*b8123daeSJed Brown 
1889d763cef2SBarry Smith .keywords: TS, timestep, solve
1890d763cef2SBarry Smith 
1891*b8123daeSJed Brown .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage()
1892d763cef2SBarry Smith @*/
1893193ac0bcSJed Brown PetscErrorCode  TSStep(TS ts)
1894d763cef2SBarry Smith {
1895362cd11cSLisandro Dalcin   PetscReal      ptime_prev;
1896dfbe8321SBarry Smith   PetscErrorCode ierr;
1897d763cef2SBarry Smith 
1898d763cef2SBarry Smith   PetscFunctionBegin;
18990700a824SBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1900d405a339SMatthew Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
1901d405a339SMatthew Knepley 
1902362cd11cSLisandro Dalcin   ts->reason = TS_CONVERGED_ITERATING;
1903362cd11cSLisandro Dalcin 
1904362cd11cSLisandro Dalcin   ptime_prev = ts->ptime;
1905d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1906193ac0bcSJed Brown   ierr = (*ts->ops->step)(ts);CHKERRQ(ierr);
1907d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr);
1908362cd11cSLisandro Dalcin   ts->time_step_prev = ts->ptime - ptime_prev;
1909362cd11cSLisandro Dalcin 
1910362cd11cSLisandro Dalcin   if (ts->reason < 0) {
1911cef5090cSJed Brown     if (ts->errorifstepfailed) {
1912cef5090cSJed Brown       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) {
1913cef5090cSJed 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]);
1914cef5090cSJed Brown       } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
1915cef5090cSJed Brown     }
1916362cd11cSLisandro Dalcin   } else if (!ts->reason) {
1917362cd11cSLisandro Dalcin     if (ts->steps >= ts->max_steps)
1918362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_ITS;
1919362cd11cSLisandro Dalcin     else if (ts->ptime >= ts->max_time)
1920362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_TIME;
1921362cd11cSLisandro Dalcin   }
1922362cd11cSLisandro Dalcin 
1923d763cef2SBarry Smith   PetscFunctionReturn(0);
1924d763cef2SBarry Smith }
1925d763cef2SBarry Smith 
19264a2ae208SSatish Balay #undef __FUNCT__
192705175c85SJed Brown #define __FUNCT__ "TSEvaluateStep"
192805175c85SJed Brown /*@
192905175c85SJed Brown    TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy.
193005175c85SJed Brown 
19311c3436cfSJed Brown    Collective on TS
193205175c85SJed Brown 
193305175c85SJed Brown    Input Arguments:
19341c3436cfSJed Brown +  ts - time stepping context
19351c3436cfSJed Brown .  order - desired order of accuracy
19361c3436cfSJed Brown -  done - whether the step was evaluated at this order (pass PETSC_NULL to generate an error if not available)
193705175c85SJed Brown 
193805175c85SJed Brown    Output Arguments:
19391c3436cfSJed Brown .  X - state at the end of the current step
194005175c85SJed Brown 
194105175c85SJed Brown    Level: advanced
194205175c85SJed Brown 
1943108c343cSJed Brown    Notes:
1944108c343cSJed Brown    This function cannot be called until all stages have been evaluated.
1945108c343cSJed 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.
1946108c343cSJed Brown 
19471c3436cfSJed Brown .seealso: TSStep(), TSAdapt
194805175c85SJed Brown @*/
19491c3436cfSJed Brown PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done)
195005175c85SJed Brown {
195105175c85SJed Brown   PetscErrorCode ierr;
195205175c85SJed Brown 
195305175c85SJed Brown   PetscFunctionBegin;
195405175c85SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
195505175c85SJed Brown   PetscValidType(ts,1);
195605175c85SJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
19571c3436cfSJed Brown   if (!ts->ops->evaluatestep) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name);
19581c3436cfSJed Brown   ierr = (*ts->ops->evaluatestep)(ts,order,X,done);CHKERRQ(ierr);
195905175c85SJed Brown   PetscFunctionReturn(0);
196005175c85SJed Brown }
196105175c85SJed Brown 
196205175c85SJed Brown #undef __FUNCT__
19636a4d4014SLisandro Dalcin #define __FUNCT__ "TSSolve"
19646a4d4014SLisandro Dalcin /*@
19656a4d4014SLisandro Dalcin    TSSolve - Steps the requested number of timesteps.
19666a4d4014SLisandro Dalcin 
19676a4d4014SLisandro Dalcin    Collective on TS
19686a4d4014SLisandro Dalcin 
19696a4d4014SLisandro Dalcin    Input Parameter:
19706a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
1971362cd11cSLisandro Dalcin -  x - the solution vector
19726a4d4014SLisandro Dalcin 
19735a3a76d0SJed Brown    Output Parameter:
19745a3a76d0SJed Brown .  ftime - time of the state vector x upon completion
19755a3a76d0SJed Brown 
19766a4d4014SLisandro Dalcin    Level: beginner
19776a4d4014SLisandro Dalcin 
19785a3a76d0SJed Brown    Notes:
19795a3a76d0SJed Brown    The final time returned by this function may be different from the time of the internally
19805a3a76d0SJed Brown    held state accessible by TSGetSolution() and TSGetTime() because the method may have
19815a3a76d0SJed Brown    stepped over the final time.
19825a3a76d0SJed Brown 
19836a4d4014SLisandro Dalcin .keywords: TS, timestep, solve
19846a4d4014SLisandro Dalcin 
19856a4d4014SLisandro Dalcin .seealso: TSCreate(), TSSetSolution(), TSStep()
19866a4d4014SLisandro Dalcin @*/
19875a3a76d0SJed Brown PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime)
19886a4d4014SLisandro Dalcin {
19894d7d938eSLisandro Dalcin   PetscBool      flg;
19904d7d938eSLisandro Dalcin   char           filename[PETSC_MAX_PATH_LEN];
19914d7d938eSLisandro Dalcin   PetscViewer    viewer;
19926a4d4014SLisandro Dalcin   PetscErrorCode ierr;
1993f22f69f0SBarry Smith 
19946a4d4014SLisandro Dalcin   PetscFunctionBegin;
19950700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1996193ac0bcSJed Brown   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
19975a3a76d0SJed Brown   if (ts->exact_final_time) {   /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */
19985a3a76d0SJed Brown     if (!ts->vec_sol || x == ts->vec_sol) {
19995a3a76d0SJed Brown       Vec y;
20005a3a76d0SJed Brown       ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
20015a3a76d0SJed Brown       ierr = TSSetSolution(ts,y);CHKERRQ(ierr);
20025a3a76d0SJed Brown       ierr = VecDestroy(&y);CHKERRQ(ierr); /* grant ownership */
20035a3a76d0SJed Brown     }
2004362cd11cSLisandro Dalcin     ierr = VecCopy(x,ts->vec_sol);CHKERRQ(ierr);
20055a3a76d0SJed Brown   } else {
2006193ac0bcSJed Brown     ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
20075a3a76d0SJed Brown   }
2008b5d403baSSean Farley   ierr = TSSetUp(ts);CHKERRQ(ierr);
20096a4d4014SLisandro Dalcin   /* reset time step and iteration counters */
2010193ac0bcSJed Brown   ts->steps = 0;
20115ef26d82SJed Brown   ts->ksp_its = 0;
20125ef26d82SJed Brown   ts->snes_its = 0;
2013c610991cSLisandro Dalcin   ts->num_snes_failures = 0;
2014c610991cSLisandro Dalcin   ts->reject = 0;
2015193ac0bcSJed Brown   ts->reason = TS_CONVERGED_ITERATING;
2016193ac0bcSJed Brown 
2017193ac0bcSJed Brown   if (ts->ops->solve) {         /* This private interface is transitional and should be removed when all implementations are updated. */
2018193ac0bcSJed Brown     ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr);
20195a3a76d0SJed Brown     ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
2020a5b82c69SSean Farley     if (ftime) *ftime = ts->ptime;
2021193ac0bcSJed Brown   } else {
20226a4d4014SLisandro Dalcin     /* steps the requested number of timesteps. */
2023362cd11cSLisandro Dalcin     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2024362cd11cSLisandro Dalcin     if (ts->steps >= ts->max_steps)
2025362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_ITS;
2026362cd11cSLisandro Dalcin     else if (ts->ptime >= ts->max_time)
2027362cd11cSLisandro Dalcin       ts->reason = TS_CONVERGED_TIME;
2028e1a7a14fSJed Brown     while (!ts->reason) {
2029193ac0bcSJed Brown       ierr = TSStep(ts);CHKERRQ(ierr);
2030193ac0bcSJed Brown       ierr = TSPostStep(ts);CHKERRQ(ierr);
2031193ac0bcSJed Brown       ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
2032193ac0bcSJed Brown     }
2033362cd11cSLisandro Dalcin     if (ts->exact_final_time && ts->ptime > ts->max_time) {
2034a43b19c4SJed Brown       ierr = TSInterpolate(ts,ts->max_time,x);CHKERRQ(ierr);
20355a3a76d0SJed Brown       if (ftime) *ftime = ts->max_time;
20360574a7fbSJed Brown     } else {
20370574a7fbSJed Brown       ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr);
20380574a7fbSJed Brown       if (ftime) *ftime = ts->ptime;
20390574a7fbSJed Brown     }
2040193ac0bcSJed Brown   }
20414d7d938eSLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
20424d7d938eSLisandro Dalcin   if (flg && !PetscPreLoadingOn) {
20434d7d938eSLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr);
20444d7d938eSLisandro Dalcin     ierr = TSView(ts,viewer);CHKERRQ(ierr);
20454d7d938eSLisandro Dalcin     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2046193ac0bcSJed Brown   }
20476a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
20486a4d4014SLisandro Dalcin }
20496a4d4014SLisandro Dalcin 
20506a4d4014SLisandro Dalcin #undef __FUNCT__
20514a2ae208SSatish Balay #define __FUNCT__ "TSMonitor"
2052228d79bcSJed Brown /*@
2053228d79bcSJed Brown    TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet()
2054228d79bcSJed Brown 
2055228d79bcSJed Brown    Collective on TS
2056228d79bcSJed Brown 
2057228d79bcSJed Brown    Input Parameters:
2058228d79bcSJed Brown +  ts - time stepping context obtained from TSCreate()
2059228d79bcSJed Brown .  step - step number that has just completed
2060228d79bcSJed Brown .  ptime - model time of the state
2061228d79bcSJed Brown -  x - state at the current model time
2062228d79bcSJed Brown 
2063228d79bcSJed Brown    Notes:
2064228d79bcSJed Brown    TSMonitor() is typically used within the time stepping implementations.
2065228d79bcSJed Brown    Users might call this function when using the TSStep() interface instead of TSSolve().
2066228d79bcSJed Brown 
2067228d79bcSJed Brown    Level: advanced
2068228d79bcSJed Brown 
2069228d79bcSJed Brown .keywords: TS, timestep
2070228d79bcSJed Brown @*/
2071a7cc72afSBarry Smith PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x)
2072d763cef2SBarry Smith {
20736849ba73SBarry Smith   PetscErrorCode ierr;
2074a7cc72afSBarry Smith   PetscInt       i,n = ts->numbermonitors;
2075d763cef2SBarry Smith 
2076d763cef2SBarry Smith   PetscFunctionBegin;
2077d763cef2SBarry Smith   for (i=0; i<n; i++) {
2078142b95e3SSatish Balay     ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr);
2079d763cef2SBarry Smith   }
2080d763cef2SBarry Smith   PetscFunctionReturn(0);
2081d763cef2SBarry Smith }
2082d763cef2SBarry Smith 
2083d763cef2SBarry Smith /* ------------------------------------------------------------------------*/
2084d763cef2SBarry Smith 
20854a2ae208SSatish Balay #undef __FUNCT__
2086a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGCreate"
2087d763cef2SBarry Smith /*@C
2088a6570f20SBarry Smith    TSMonitorLGCreate - Creates a line graph context for use with
2089d763cef2SBarry Smith    TS to monitor convergence of preconditioned residual norms.
2090d763cef2SBarry Smith 
2091d763cef2SBarry Smith    Collective on TS
2092d763cef2SBarry Smith 
2093d763cef2SBarry Smith    Input Parameters:
2094d763cef2SBarry Smith +  host - the X display to open, or null for the local machine
2095d763cef2SBarry Smith .  label - the title to put in the title bar
20967c922b88SBarry Smith .  x, y - the screen coordinates of the upper left coordinate of the window
2097d763cef2SBarry Smith -  m, n - the screen width and height in pixels
2098d763cef2SBarry Smith 
2099d763cef2SBarry Smith    Output Parameter:
2100d763cef2SBarry Smith .  draw - the drawing context
2101d763cef2SBarry Smith 
2102d763cef2SBarry Smith    Options Database Key:
2103a6570f20SBarry Smith .  -ts_monitor_draw - automatically sets line graph monitor
2104d763cef2SBarry Smith 
2105d763cef2SBarry Smith    Notes:
2106a6570f20SBarry Smith    Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy().
2107d763cef2SBarry Smith 
2108d763cef2SBarry Smith    Level: intermediate
2109d763cef2SBarry Smith 
21107c922b88SBarry Smith .keywords: TS, monitor, line graph, residual, seealso
2111d763cef2SBarry Smith 
2112a6570f20SBarry Smith .seealso: TSMonitorLGDestroy(), TSMonitorSet()
21137c922b88SBarry Smith 
2114d763cef2SBarry Smith @*/
21157087cfbeSBarry Smith PetscErrorCode  TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2116d763cef2SBarry Smith {
2117b0a32e0cSBarry Smith   PetscDraw      win;
2118dfbe8321SBarry Smith   PetscErrorCode ierr;
2119d763cef2SBarry Smith 
2120d763cef2SBarry Smith   PetscFunctionBegin;
2121b0a32e0cSBarry Smith   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
2122b0a32e0cSBarry Smith   ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr);
2123b0a32e0cSBarry Smith   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
2124b0a32e0cSBarry Smith   ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr);
2125d763cef2SBarry Smith 
212652e6d16bSBarry Smith   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
2127d763cef2SBarry Smith   PetscFunctionReturn(0);
2128d763cef2SBarry Smith }
2129d763cef2SBarry Smith 
21304a2ae208SSatish Balay #undef __FUNCT__
2131a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLG"
2132a6570f20SBarry Smith PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx)
2133d763cef2SBarry Smith {
2134b0a32e0cSBarry Smith   PetscDrawLG    lg = (PetscDrawLG) monctx;
213587828ca2SBarry Smith   PetscReal      x,y = ptime;
2136dfbe8321SBarry Smith   PetscErrorCode ierr;
2137d763cef2SBarry Smith 
2138d763cef2SBarry Smith   PetscFunctionBegin;
21397c922b88SBarry Smith   if (!monctx) {
21407c922b88SBarry Smith     MPI_Comm    comm;
2141b0a32e0cSBarry Smith     PetscViewer viewer;
21427c922b88SBarry Smith 
21437c922b88SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
2144b0a32e0cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(comm);
2145b0a32e0cSBarry Smith     ierr   = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
21467c922b88SBarry Smith   }
21477c922b88SBarry Smith 
2148b0a32e0cSBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
214987828ca2SBarry Smith   x = (PetscReal)n;
2150b0a32e0cSBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2151d763cef2SBarry Smith   if (n < 20 || (n % 5)) {
2152b0a32e0cSBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2153d763cef2SBarry Smith   }
2154d763cef2SBarry Smith   PetscFunctionReturn(0);
2155d763cef2SBarry Smith }
2156d763cef2SBarry Smith 
21574a2ae208SSatish Balay #undef __FUNCT__
2158a6570f20SBarry Smith #define __FUNCT__ "TSMonitorLGDestroy"
2159d763cef2SBarry Smith /*@C
2160a6570f20SBarry Smith    TSMonitorLGDestroy - Destroys a line graph context that was created
2161a6570f20SBarry Smith    with TSMonitorLGCreate().
2162d763cef2SBarry Smith 
2163b0a32e0cSBarry Smith    Collective on PetscDrawLG
2164d763cef2SBarry Smith 
2165d763cef2SBarry Smith    Input Parameter:
2166d763cef2SBarry Smith .  draw - the drawing context
2167d763cef2SBarry Smith 
2168d763cef2SBarry Smith    Level: intermediate
2169d763cef2SBarry Smith 
2170d763cef2SBarry Smith .keywords: TS, monitor, line graph, destroy
2171d763cef2SBarry Smith 
2172a6570f20SBarry Smith .seealso: TSMonitorLGCreate(),  TSMonitorSet(), TSMonitorLG();
2173d763cef2SBarry Smith @*/
21746bf464f9SBarry Smith PetscErrorCode  TSMonitorLGDestroy(PetscDrawLG *drawlg)
2175d763cef2SBarry Smith {
2176b0a32e0cSBarry Smith   PetscDraw      draw;
2177dfbe8321SBarry Smith   PetscErrorCode ierr;
2178d763cef2SBarry Smith 
2179d763cef2SBarry Smith   PetscFunctionBegin;
21806bf464f9SBarry Smith   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
21816bf464f9SBarry Smith   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
2182b0a32e0cSBarry Smith   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
2183d763cef2SBarry Smith   PetscFunctionReturn(0);
2184d763cef2SBarry Smith }
2185d763cef2SBarry Smith 
21864a2ae208SSatish Balay #undef __FUNCT__
21874a2ae208SSatish Balay #define __FUNCT__ "TSGetTime"
2188d763cef2SBarry Smith /*@
2189*b8123daeSJed Brown    TSGetTime - Gets the time of the most recently completed step.
2190d763cef2SBarry Smith 
2191d763cef2SBarry Smith    Not Collective
2192d763cef2SBarry Smith 
2193d763cef2SBarry Smith    Input Parameter:
2194d763cef2SBarry Smith .  ts - the TS context obtained from TSCreate()
2195d763cef2SBarry Smith 
2196d763cef2SBarry Smith    Output Parameter:
2197d763cef2SBarry Smith .  t  - the current time
2198d763cef2SBarry Smith 
2199d763cef2SBarry Smith    Level: beginner
2200d763cef2SBarry Smith 
2201*b8123daeSJed Brown    Note:
2202*b8123daeSJed Brown    When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(),
2203*b8123daeSJed Brown    TSSetPreStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated.
2204*b8123daeSJed Brown 
2205d763cef2SBarry Smith .seealso: TSSetInitialTimeStep(), TSGetTimeStep()
2206d763cef2SBarry Smith 
2207d763cef2SBarry Smith .keywords: TS, get, time
2208d763cef2SBarry Smith @*/
22097087cfbeSBarry Smith PetscErrorCode  TSGetTime(TS ts,PetscReal* t)
2210d763cef2SBarry Smith {
2211d763cef2SBarry Smith   PetscFunctionBegin;
22120700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
22134482741eSBarry Smith   PetscValidDoublePointer(t,2);
2214d763cef2SBarry Smith   *t = ts->ptime;
2215d763cef2SBarry Smith   PetscFunctionReturn(0);
2216d763cef2SBarry Smith }
2217d763cef2SBarry Smith 
22184a2ae208SSatish Balay #undef __FUNCT__
22196a4d4014SLisandro Dalcin #define __FUNCT__ "TSSetTime"
22206a4d4014SLisandro Dalcin /*@
22216a4d4014SLisandro Dalcin    TSSetTime - Allows one to reset the time.
22226a4d4014SLisandro Dalcin 
22233f9fe445SBarry Smith    Logically Collective on TS
22246a4d4014SLisandro Dalcin 
22256a4d4014SLisandro Dalcin    Input Parameters:
22266a4d4014SLisandro Dalcin +  ts - the TS context obtained from TSCreate()
22276a4d4014SLisandro Dalcin -  time - the time
22286a4d4014SLisandro Dalcin 
22296a4d4014SLisandro Dalcin    Level: intermediate
22306a4d4014SLisandro Dalcin 
22316a4d4014SLisandro Dalcin .seealso: TSGetTime(), TSSetDuration()
22326a4d4014SLisandro Dalcin 
22336a4d4014SLisandro Dalcin .keywords: TS, set, time
22346a4d4014SLisandro Dalcin @*/
22357087cfbeSBarry Smith PetscErrorCode  TSSetTime(TS ts, PetscReal t)
22366a4d4014SLisandro Dalcin {
22376a4d4014SLisandro Dalcin   PetscFunctionBegin;
22380700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2239c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,t,2);
22406a4d4014SLisandro Dalcin   ts->ptime = t;
22416a4d4014SLisandro Dalcin   PetscFunctionReturn(0);
22426a4d4014SLisandro Dalcin }
22436a4d4014SLisandro Dalcin 
22446a4d4014SLisandro Dalcin #undef __FUNCT__
22454a2ae208SSatish Balay #define __FUNCT__ "TSSetOptionsPrefix"
2246d763cef2SBarry Smith /*@C
2247d763cef2SBarry Smith    TSSetOptionsPrefix - Sets the prefix used for searching for all
2248d763cef2SBarry Smith    TS options in the database.
2249d763cef2SBarry Smith 
22503f9fe445SBarry Smith    Logically Collective on TS
2251d763cef2SBarry Smith 
2252d763cef2SBarry Smith    Input Parameter:
2253d763cef2SBarry Smith +  ts     - The TS context
2254d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
2255d763cef2SBarry Smith 
2256d763cef2SBarry Smith    Notes:
2257d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
2258d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
2259d763cef2SBarry Smith    hyphen.
2260d763cef2SBarry Smith 
2261d763cef2SBarry Smith    Level: advanced
2262d763cef2SBarry Smith 
2263d763cef2SBarry Smith .keywords: TS, set, options, prefix, database
2264d763cef2SBarry Smith 
2265d763cef2SBarry Smith .seealso: TSSetFromOptions()
2266d763cef2SBarry Smith 
2267d763cef2SBarry Smith @*/
22687087cfbeSBarry Smith PetscErrorCode  TSSetOptionsPrefix(TS ts,const char prefix[])
2269d763cef2SBarry Smith {
2270dfbe8321SBarry Smith   PetscErrorCode ierr;
2271089b2837SJed Brown   SNES           snes;
2272d763cef2SBarry Smith 
2273d763cef2SBarry Smith   PetscFunctionBegin;
22740700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2275d763cef2SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2276089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2277089b2837SJed Brown   ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2278d763cef2SBarry Smith   PetscFunctionReturn(0);
2279d763cef2SBarry Smith }
2280d763cef2SBarry Smith 
2281d763cef2SBarry Smith 
22824a2ae208SSatish Balay #undef __FUNCT__
22834a2ae208SSatish Balay #define __FUNCT__ "TSAppendOptionsPrefix"
2284d763cef2SBarry Smith /*@C
2285d763cef2SBarry Smith    TSAppendOptionsPrefix - Appends to the prefix used for searching for all
2286d763cef2SBarry Smith    TS options in the database.
2287d763cef2SBarry Smith 
22883f9fe445SBarry Smith    Logically Collective on TS
2289d763cef2SBarry Smith 
2290d763cef2SBarry Smith    Input Parameter:
2291d763cef2SBarry Smith +  ts     - The TS context
2292d763cef2SBarry Smith -  prefix - The prefix to prepend to all option names
2293d763cef2SBarry Smith 
2294d763cef2SBarry Smith    Notes:
2295d763cef2SBarry Smith    A hyphen (-) must NOT be given at the beginning of the prefix name.
2296d763cef2SBarry Smith    The first character of all runtime options is AUTOMATICALLY the
2297d763cef2SBarry Smith    hyphen.
2298d763cef2SBarry Smith 
2299d763cef2SBarry Smith    Level: advanced
2300d763cef2SBarry Smith 
2301d763cef2SBarry Smith .keywords: TS, append, options, prefix, database
2302d763cef2SBarry Smith 
2303d763cef2SBarry Smith .seealso: TSGetOptionsPrefix()
2304d763cef2SBarry Smith 
2305d763cef2SBarry Smith @*/
23067087cfbeSBarry Smith PetscErrorCode  TSAppendOptionsPrefix(TS ts,const char prefix[])
2307d763cef2SBarry Smith {
2308dfbe8321SBarry Smith   PetscErrorCode ierr;
2309089b2837SJed Brown   SNES           snes;
2310d763cef2SBarry Smith 
2311d763cef2SBarry Smith   PetscFunctionBegin;
23120700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2313d763cef2SBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2314089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2315089b2837SJed Brown   ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr);
2316d763cef2SBarry Smith   PetscFunctionReturn(0);
2317d763cef2SBarry Smith }
2318d763cef2SBarry Smith 
23194a2ae208SSatish Balay #undef __FUNCT__
23204a2ae208SSatish Balay #define __FUNCT__ "TSGetOptionsPrefix"
2321d763cef2SBarry Smith /*@C
2322d763cef2SBarry Smith    TSGetOptionsPrefix - Sets the prefix used for searching for all
2323d763cef2SBarry Smith    TS options in the database.
2324d763cef2SBarry Smith 
2325d763cef2SBarry Smith    Not Collective
2326d763cef2SBarry Smith 
2327d763cef2SBarry Smith    Input Parameter:
2328d763cef2SBarry Smith .  ts - The TS context
2329d763cef2SBarry Smith 
2330d763cef2SBarry Smith    Output Parameter:
2331d763cef2SBarry Smith .  prefix - A pointer to the prefix string used
2332d763cef2SBarry Smith 
2333d763cef2SBarry Smith    Notes: On the fortran side, the user should pass in a string 'prifix' of
2334d763cef2SBarry Smith    sufficient length to hold the prefix.
2335d763cef2SBarry Smith 
2336d763cef2SBarry Smith    Level: intermediate
2337d763cef2SBarry Smith 
2338d763cef2SBarry Smith .keywords: TS, get, options, prefix, database
2339d763cef2SBarry Smith 
2340d763cef2SBarry Smith .seealso: TSAppendOptionsPrefix()
2341d763cef2SBarry Smith @*/
23427087cfbeSBarry Smith PetscErrorCode  TSGetOptionsPrefix(TS ts,const char *prefix[])
2343d763cef2SBarry Smith {
2344dfbe8321SBarry Smith   PetscErrorCode ierr;
2345d763cef2SBarry Smith 
2346d763cef2SBarry Smith   PetscFunctionBegin;
23470700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
23484482741eSBarry Smith   PetscValidPointer(prefix,2);
2349d763cef2SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr);
2350d763cef2SBarry Smith   PetscFunctionReturn(0);
2351d763cef2SBarry Smith }
2352d763cef2SBarry Smith 
23534a2ae208SSatish Balay #undef __FUNCT__
23544a2ae208SSatish Balay #define __FUNCT__ "TSGetRHSJacobian"
2355d763cef2SBarry Smith /*@C
2356d763cef2SBarry Smith    TSGetRHSJacobian - Returns the Jacobian J at the present timestep.
2357d763cef2SBarry Smith 
2358d763cef2SBarry Smith    Not Collective, but parallel objects are returned if TS is parallel
2359d763cef2SBarry Smith 
2360d763cef2SBarry Smith    Input Parameter:
2361d763cef2SBarry Smith .  ts  - The TS context obtained from TSCreate()
2362d763cef2SBarry Smith 
2363d763cef2SBarry Smith    Output Parameters:
2364d763cef2SBarry Smith +  J   - The Jacobian J of F, where U_t = F(U,t)
2365d763cef2SBarry Smith .  M   - The preconditioner matrix, usually the same as J
2366089b2837SJed Brown .  func - Function to compute the Jacobian of the RHS
2367d763cef2SBarry Smith -  ctx - User-defined context for Jacobian evaluation routine
2368d763cef2SBarry Smith 
2369d763cef2SBarry Smith    Notes: You can pass in PETSC_NULL for any return argument you do not need.
2370d763cef2SBarry Smith 
2371d763cef2SBarry Smith    Level: intermediate
2372d763cef2SBarry Smith 
237326d46c62SHong Zhang .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
2374d763cef2SBarry Smith 
2375d763cef2SBarry Smith .keywords: TS, timestep, get, matrix, Jacobian
2376d763cef2SBarry Smith @*/
2377089b2837SJed Brown PetscErrorCode  TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx)
2378d763cef2SBarry Smith {
2379089b2837SJed Brown   PetscErrorCode ierr;
2380089b2837SJed Brown   SNES           snes;
2381089b2837SJed Brown 
2382d763cef2SBarry Smith   PetscFunctionBegin;
2383089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2384089b2837SJed Brown   ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
23853b5f76d0SSean Farley   if (func) *func = ts->userops->rhsjacobian;
238626d46c62SHong Zhang   if (ctx) *ctx = ts->jacP;
2387d763cef2SBarry Smith   PetscFunctionReturn(0);
2388d763cef2SBarry Smith }
2389d763cef2SBarry Smith 
23901713a123SBarry Smith #undef __FUNCT__
23912eca1d9cSJed Brown #define __FUNCT__ "TSGetIJacobian"
23922eca1d9cSJed Brown /*@C
23932eca1d9cSJed Brown    TSGetIJacobian - Returns the implicit Jacobian at the present timestep.
23942eca1d9cSJed Brown 
23952eca1d9cSJed Brown    Not Collective, but parallel objects are returned if TS is parallel
23962eca1d9cSJed Brown 
23972eca1d9cSJed Brown    Input Parameter:
23982eca1d9cSJed Brown .  ts  - The TS context obtained from TSCreate()
23992eca1d9cSJed Brown 
24002eca1d9cSJed Brown    Output Parameters:
24012eca1d9cSJed Brown +  A   - The Jacobian of F(t,U,U_t)
24022eca1d9cSJed Brown .  B   - The preconditioner matrix, often the same as A
24032eca1d9cSJed Brown .  f   - The function to compute the matrices
24042eca1d9cSJed Brown - ctx - User-defined context for Jacobian evaluation routine
24052eca1d9cSJed Brown 
24062eca1d9cSJed Brown    Notes: You can pass in PETSC_NULL for any return argument you do not need.
24072eca1d9cSJed Brown 
24082eca1d9cSJed Brown    Level: advanced
24092eca1d9cSJed Brown 
24102eca1d9cSJed Brown .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber()
24112eca1d9cSJed Brown 
24122eca1d9cSJed Brown .keywords: TS, timestep, get, matrix, Jacobian
24132eca1d9cSJed Brown @*/
24147087cfbeSBarry Smith PetscErrorCode  TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx)
24152eca1d9cSJed Brown {
2416089b2837SJed Brown   PetscErrorCode ierr;
2417089b2837SJed Brown   SNES           snes;
2418089b2837SJed Brown 
24192eca1d9cSJed Brown   PetscFunctionBegin;
2420089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2421089b2837SJed Brown   ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
24223b5f76d0SSean Farley   if (f) *f = ts->userops->ijacobian;
24232eca1d9cSJed Brown   if (ctx) *ctx = ts->jacP;
24242eca1d9cSJed Brown   PetscFunctionReturn(0);
24252eca1d9cSJed Brown }
24262eca1d9cSJed Brown 
24276083293cSBarry Smith typedef struct {
24286083293cSBarry Smith   PetscViewer viewer;
24296083293cSBarry Smith   Vec         initialsolution;
24306083293cSBarry Smith   PetscBool   showinitial;
24316083293cSBarry Smith } TSMonitorSolutionCtx;
24326083293cSBarry Smith 
24332eca1d9cSJed Brown #undef __FUNCT__
2434a6570f20SBarry Smith #define __FUNCT__ "TSMonitorSolution"
24351713a123SBarry Smith /*@C
2436a6570f20SBarry Smith    TSMonitorSolution - Monitors progress of the TS solvers by calling
24371713a123SBarry Smith    VecView() for the solution at each timestep
24381713a123SBarry Smith 
24391713a123SBarry Smith    Collective on TS
24401713a123SBarry Smith 
24411713a123SBarry Smith    Input Parameters:
24421713a123SBarry Smith +  ts - the TS context
24431713a123SBarry Smith .  step - current time-step
2444142b95e3SSatish Balay .  ptime - current time
24451713a123SBarry Smith -  dummy - either a viewer or PETSC_NULL
24461713a123SBarry Smith 
24471713a123SBarry Smith    Level: intermediate
24481713a123SBarry Smith 
24491713a123SBarry Smith .keywords: TS,  vector, monitor, view
24501713a123SBarry Smith 
2451a6570f20SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
24521713a123SBarry Smith @*/
24537087cfbeSBarry Smith PetscErrorCode  TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy)
24541713a123SBarry Smith {
2455dfbe8321SBarry Smith   PetscErrorCode       ierr;
24566083293cSBarry Smith   TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy;
24571713a123SBarry Smith 
24581713a123SBarry Smith   PetscFunctionBegin;
24596083293cSBarry Smith   if (!step && ictx->showinitial) {
24606083293cSBarry Smith     if (!ictx->initialsolution) {
24616083293cSBarry Smith       ierr = VecDuplicate(x,&ictx->initialsolution);CHKERRQ(ierr);
24621713a123SBarry Smith     }
24636083293cSBarry Smith     ierr = VecCopy(x,ictx->initialsolution);CHKERRQ(ierr);
24646083293cSBarry Smith   }
24656083293cSBarry Smith   if (ictx->showinitial) {
24666083293cSBarry Smith     PetscReal pause;
24676083293cSBarry Smith     ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr);
24686083293cSBarry Smith     ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr);
24696083293cSBarry Smith     ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr);
24706083293cSBarry Smith     ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr);
24716083293cSBarry Smith     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr);
24726083293cSBarry Smith   }
24736083293cSBarry Smith   ierr = VecView(x,ictx->viewer);CHKERRQ(ierr);
24746083293cSBarry Smith   if (ictx->showinitial) {
24756083293cSBarry Smith     ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr);
24766083293cSBarry Smith   }
24771713a123SBarry Smith   PetscFunctionReturn(0);
24781713a123SBarry Smith }
24791713a123SBarry Smith 
24801713a123SBarry Smith 
24816c699258SBarry Smith #undef __FUNCT__
24826083293cSBarry Smith #define __FUNCT__ "TSMonitorSolutionDestroy"
24836083293cSBarry Smith /*@C
24846083293cSBarry Smith    TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution
24856083293cSBarry Smith 
24866083293cSBarry Smith    Collective on TS
24876083293cSBarry Smith 
24886083293cSBarry Smith    Input Parameters:
24896083293cSBarry Smith .    ctx - the monitor context
24906083293cSBarry Smith 
24916083293cSBarry Smith    Level: intermediate
24926083293cSBarry Smith 
24936083293cSBarry Smith .keywords: TS,  vector, monitor, view
24946083293cSBarry Smith 
24956083293cSBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
24966083293cSBarry Smith @*/
24976083293cSBarry Smith PetscErrorCode  TSMonitorSolutionDestroy(void **ctx)
24986083293cSBarry Smith {
24996083293cSBarry Smith   PetscErrorCode       ierr;
25006083293cSBarry Smith   TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx;
25016083293cSBarry Smith 
25026083293cSBarry Smith   PetscFunctionBegin;
25036083293cSBarry Smith   ierr = PetscViewerDestroy(&ictx->viewer);CHKERRQ(ierr);
25046083293cSBarry Smith   ierr = VecDestroy(&ictx->initialsolution);CHKERRQ(ierr);
25056083293cSBarry Smith   ierr = PetscFree(ictx);CHKERRQ(ierr);
25066083293cSBarry Smith   PetscFunctionReturn(0);
25076083293cSBarry Smith }
25086083293cSBarry Smith 
25096083293cSBarry Smith #undef __FUNCT__
25106083293cSBarry Smith #define __FUNCT__ "TSMonitorSolutionCreate"
25116083293cSBarry Smith /*@C
25126083293cSBarry Smith    TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution
25136083293cSBarry Smith 
25146083293cSBarry Smith    Collective on TS
25156083293cSBarry Smith 
25166083293cSBarry Smith    Input Parameter:
25176083293cSBarry Smith .    ts - time-step context
25186083293cSBarry Smith 
25196083293cSBarry Smith    Output Patameter:
25206083293cSBarry Smith .    ctx - the monitor context
25216083293cSBarry Smith 
25226083293cSBarry Smith    Level: intermediate
25236083293cSBarry Smith 
25246083293cSBarry Smith .keywords: TS,  vector, monitor, view
25256083293cSBarry Smith 
25266083293cSBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution()
25276083293cSBarry Smith @*/
25286083293cSBarry Smith PetscErrorCode  TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx)
25296083293cSBarry Smith {
25306083293cSBarry Smith   PetscErrorCode       ierr;
25316083293cSBarry Smith   TSMonitorSolutionCtx *ictx;
25326083293cSBarry Smith 
25336083293cSBarry Smith   PetscFunctionBegin;
25346083293cSBarry Smith   ierr = PetscNew(TSMonitorSolutionCtx,&ictx);CHKERRQ(ierr);
25356083293cSBarry Smith   *ctx = (void*)ictx;
25366083293cSBarry Smith   if (!viewer) {
25376083293cSBarry Smith     viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm);
25386083293cSBarry Smith   }
25396083293cSBarry Smith   ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
25406083293cSBarry Smith   ictx->viewer      = viewer;
25416083293cSBarry Smith   ictx->showinitial = PETSC_FALSE;
25426083293cSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);CHKERRQ(ierr);
25436083293cSBarry Smith   PetscFunctionReturn(0);
25446083293cSBarry Smith }
25456083293cSBarry Smith 
25466083293cSBarry Smith #undef __FUNCT__
25476c699258SBarry Smith #define __FUNCT__ "TSSetDM"
25486c699258SBarry Smith /*@
25496c699258SBarry Smith    TSSetDM - Sets the DM that may be used by some preconditioners
25506c699258SBarry Smith 
25513f9fe445SBarry Smith    Logically Collective on TS and DM
25526c699258SBarry Smith 
25536c699258SBarry Smith    Input Parameters:
25546c699258SBarry Smith +  ts - the preconditioner context
25556c699258SBarry Smith -  dm - the dm
25566c699258SBarry Smith 
25576c699258SBarry Smith    Level: intermediate
25586c699258SBarry Smith 
25596c699258SBarry Smith 
25606c699258SBarry Smith .seealso: TSGetDM(), SNESSetDM(), SNESGetDM()
25616c699258SBarry Smith @*/
25627087cfbeSBarry Smith PetscErrorCode  TSSetDM(TS ts,DM dm)
25636c699258SBarry Smith {
25646c699258SBarry Smith   PetscErrorCode ierr;
2565089b2837SJed Brown   SNES           snes;
25666c699258SBarry Smith 
25676c699258SBarry Smith   PetscFunctionBegin;
25680700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
256970663e4aSLisandro Dalcin   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
25706bf464f9SBarry Smith   ierr = DMDestroy(&ts->dm);CHKERRQ(ierr);
25716c699258SBarry Smith   ts->dm = dm;
2572089b2837SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2573089b2837SJed Brown   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
25746c699258SBarry Smith   PetscFunctionReturn(0);
25756c699258SBarry Smith }
25766c699258SBarry Smith 
25776c699258SBarry Smith #undef __FUNCT__
25786c699258SBarry Smith #define __FUNCT__ "TSGetDM"
25796c699258SBarry Smith /*@
25806c699258SBarry Smith    TSGetDM - Gets the DM that may be used by some preconditioners
25816c699258SBarry Smith 
25823f9fe445SBarry Smith    Not Collective
25836c699258SBarry Smith 
25846c699258SBarry Smith    Input Parameter:
25856c699258SBarry Smith . ts - the preconditioner context
25866c699258SBarry Smith 
25876c699258SBarry Smith    Output Parameter:
25886c699258SBarry Smith .  dm - the dm
25896c699258SBarry Smith 
25906c699258SBarry Smith    Level: intermediate
25916c699258SBarry Smith 
25926c699258SBarry Smith 
25936c699258SBarry Smith .seealso: TSSetDM(), SNESSetDM(), SNESGetDM()
25946c699258SBarry Smith @*/
25957087cfbeSBarry Smith PetscErrorCode  TSGetDM(TS ts,DM *dm)
25966c699258SBarry Smith {
2597496e6a7aSJed Brown   PetscErrorCode ierr;
2598496e6a7aSJed Brown 
25996c699258SBarry Smith   PetscFunctionBegin;
26000700a824SBarry Smith   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2601496e6a7aSJed Brown   if (!ts->dm) {
2602496e6a7aSJed Brown     ierr = DMShellCreate(((PetscObject)ts)->comm,&ts->dm);CHKERRQ(ierr);
2603496e6a7aSJed Brown     if (ts->snes) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);}
2604496e6a7aSJed Brown   }
26056c699258SBarry Smith   *dm = ts->dm;
26066c699258SBarry Smith   PetscFunctionReturn(0);
26076c699258SBarry Smith }
26081713a123SBarry Smith 
26090f5c6efeSJed Brown #undef __FUNCT__
26100f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction"
26110f5c6efeSJed Brown /*@
26120f5c6efeSJed Brown    SNESTSFormFunction - Function to evaluate nonlinear residual
26130f5c6efeSJed Brown 
26143f9fe445SBarry Smith    Logically Collective on SNES
26150f5c6efeSJed Brown 
26160f5c6efeSJed Brown    Input Parameter:
2617d42a1c89SJed Brown + snes - nonlinear solver
26180f5c6efeSJed Brown . X - the current state at which to evaluate the residual
2619d42a1c89SJed Brown - ctx - user context, must be a TS
26200f5c6efeSJed Brown 
26210f5c6efeSJed Brown    Output Parameter:
26220f5c6efeSJed Brown . F - the nonlinear residual
26230f5c6efeSJed Brown 
26240f5c6efeSJed Brown    Notes:
26250f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
26260f5c6efeSJed Brown    It is most frequently passed to MatFDColoringSetFunction().
26270f5c6efeSJed Brown 
26280f5c6efeSJed Brown    Level: advanced
26290f5c6efeSJed Brown 
26300f5c6efeSJed Brown .seealso: SNESSetFunction(), MatFDColoringSetFunction()
26310f5c6efeSJed Brown @*/
26327087cfbeSBarry Smith PetscErrorCode  SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx)
26330f5c6efeSJed Brown {
26340f5c6efeSJed Brown   TS ts = (TS)ctx;
26350f5c6efeSJed Brown   PetscErrorCode ierr;
26360f5c6efeSJed Brown 
26370f5c6efeSJed Brown   PetscFunctionBegin;
26380f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
26390f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
26400f5c6efeSJed Brown   PetscValidHeaderSpecific(F,VEC_CLASSID,3);
26410f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,4);
26420f5c6efeSJed Brown   ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr);
26430f5c6efeSJed Brown   PetscFunctionReturn(0);
26440f5c6efeSJed Brown }
26450f5c6efeSJed Brown 
26460f5c6efeSJed Brown #undef __FUNCT__
26470f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian"
26480f5c6efeSJed Brown /*@
26490f5c6efeSJed Brown    SNESTSFormJacobian - Function to evaluate the Jacobian
26500f5c6efeSJed Brown 
26510f5c6efeSJed Brown    Collective on SNES
26520f5c6efeSJed Brown 
26530f5c6efeSJed Brown    Input Parameter:
26540f5c6efeSJed Brown + snes - nonlinear solver
26550f5c6efeSJed Brown . X - the current state at which to evaluate the residual
26560f5c6efeSJed Brown - ctx - user context, must be a TS
26570f5c6efeSJed Brown 
26580f5c6efeSJed Brown    Output Parameter:
26590f5c6efeSJed Brown + A - the Jacobian
26600f5c6efeSJed Brown . B - the preconditioning matrix (may be the same as A)
26610f5c6efeSJed Brown - flag - indicates any structure change in the matrix
26620f5c6efeSJed Brown 
26630f5c6efeSJed Brown    Notes:
26640f5c6efeSJed Brown    This function is not normally called by users and is automatically registered with the SNES used by TS.
26650f5c6efeSJed Brown 
26660f5c6efeSJed Brown    Level: developer
26670f5c6efeSJed Brown 
26680f5c6efeSJed Brown .seealso: SNESSetJacobian()
26690f5c6efeSJed Brown @*/
26707087cfbeSBarry Smith PetscErrorCode  SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx)
26710f5c6efeSJed Brown {
26720f5c6efeSJed Brown   TS ts = (TS)ctx;
26730f5c6efeSJed Brown   PetscErrorCode ierr;
26740f5c6efeSJed Brown 
26750f5c6efeSJed Brown   PetscFunctionBegin;
26760f5c6efeSJed Brown   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
26770f5c6efeSJed Brown   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
26780f5c6efeSJed Brown   PetscValidPointer(A,3);
26790f5c6efeSJed Brown   PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
26800f5c6efeSJed Brown   PetscValidPointer(B,4);
26810f5c6efeSJed Brown   PetscValidHeaderSpecific(*B,MAT_CLASSID,4);
26820f5c6efeSJed Brown   PetscValidPointer(flag,5);
26830f5c6efeSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,6);
26840f5c6efeSJed Brown   ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr);
26850f5c6efeSJed Brown   PetscFunctionReturn(0);
26860f5c6efeSJed Brown }
2687325fc9f4SBarry Smith 
26880e4ef248SJed Brown #undef __FUNCT__
26890e4ef248SJed Brown #define __FUNCT__ "TSComputeRHSFunctionLinear"
26900e4ef248SJed Brown /*@C
26910e4ef248SJed Brown    TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only
26920e4ef248SJed Brown 
26930e4ef248SJed Brown    Collective on TS
26940e4ef248SJed Brown 
26950e4ef248SJed Brown    Input Arguments:
26960e4ef248SJed Brown +  ts - time stepping context
26970e4ef248SJed Brown .  t - time at which to evaluate
26980e4ef248SJed Brown .  X - state at which to evaluate
26990e4ef248SJed Brown -  ctx - context
27000e4ef248SJed Brown 
27010e4ef248SJed Brown    Output Arguments:
27020e4ef248SJed Brown .  F - right hand side
27030e4ef248SJed Brown 
27040e4ef248SJed Brown    Level: intermediate
27050e4ef248SJed Brown 
27060e4ef248SJed Brown    Notes:
27070e4ef248SJed Brown    This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems.
27080e4ef248SJed Brown    The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian().
27090e4ef248SJed Brown 
27100e4ef248SJed Brown .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant()
27110e4ef248SJed Brown @*/
27120e4ef248SJed Brown PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
27130e4ef248SJed Brown {
27140e4ef248SJed Brown   PetscErrorCode ierr;
27150e4ef248SJed Brown   Mat Arhs,Brhs;
27160e4ef248SJed Brown   MatStructure flg2;
27170e4ef248SJed Brown 
27180e4ef248SJed Brown   PetscFunctionBegin;
27190e4ef248SJed Brown   ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr);
27200e4ef248SJed Brown   ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr);
27210e4ef248SJed Brown   ierr = MatMult(Arhs,X,F);CHKERRQ(ierr);
27220e4ef248SJed Brown   PetscFunctionReturn(0);
27230e4ef248SJed Brown }
27240e4ef248SJed Brown 
27250e4ef248SJed Brown #undef __FUNCT__
27260e4ef248SJed Brown #define __FUNCT__ "TSComputeRHSJacobianConstant"
27270e4ef248SJed Brown /*@C
27280e4ef248SJed Brown    TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent.
27290e4ef248SJed Brown 
27300e4ef248SJed Brown    Collective on TS
27310e4ef248SJed Brown 
27320e4ef248SJed Brown    Input Arguments:
27330e4ef248SJed Brown +  ts - time stepping context
27340e4ef248SJed Brown .  t - time at which to evaluate
27350e4ef248SJed Brown .  X - state at which to evaluate
27360e4ef248SJed Brown -  ctx - context
27370e4ef248SJed Brown 
27380e4ef248SJed Brown    Output Arguments:
27390e4ef248SJed Brown +  A - pointer to operator
27400e4ef248SJed Brown .  B - pointer to preconditioning matrix
27410e4ef248SJed Brown -  flg - matrix structure flag
27420e4ef248SJed Brown 
27430e4ef248SJed Brown    Level: intermediate
27440e4ef248SJed Brown 
27450e4ef248SJed Brown    Notes:
27460e4ef248SJed Brown    This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems.
27470e4ef248SJed Brown 
27480e4ef248SJed Brown .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear()
27490e4ef248SJed Brown @*/
27500e4ef248SJed Brown PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx)
27510e4ef248SJed Brown {
27520e4ef248SJed Brown 
27530e4ef248SJed Brown   PetscFunctionBegin;
27540e4ef248SJed Brown   *flg = SAME_PRECONDITIONER;
27550e4ef248SJed Brown   PetscFunctionReturn(0);
27560e4ef248SJed Brown }
27570e4ef248SJed Brown 
27580026cea9SSean Farley #undef __FUNCT__
27590026cea9SSean Farley #define __FUNCT__ "TSComputeIFunctionLinear"
27600026cea9SSean Farley /*@C
27610026cea9SSean Farley    TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only
27620026cea9SSean Farley 
27630026cea9SSean Farley    Collective on TS
27640026cea9SSean Farley 
27650026cea9SSean Farley    Input Arguments:
27660026cea9SSean Farley +  ts - time stepping context
27670026cea9SSean Farley .  t - time at which to evaluate
27680026cea9SSean Farley .  X - state at which to evaluate
27690026cea9SSean Farley .  Xdot - time derivative of state vector
27700026cea9SSean Farley -  ctx - context
27710026cea9SSean Farley 
27720026cea9SSean Farley    Output Arguments:
27730026cea9SSean Farley .  F - left hand side
27740026cea9SSean Farley 
27750026cea9SSean Farley    Level: intermediate
27760026cea9SSean Farley 
27770026cea9SSean Farley    Notes:
27780026cea9SSean 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
27790026cea9SSean Farley    user is required to write their own TSComputeIFunction.
27800026cea9SSean Farley    This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems.
27810026cea9SSean Farley    The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian().
27820026cea9SSean Farley 
27830026cea9SSean Farley .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant()
27840026cea9SSean Farley @*/
27850026cea9SSean Farley PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
27860026cea9SSean Farley {
27870026cea9SSean Farley   PetscErrorCode ierr;
27880026cea9SSean Farley   Mat A,B;
27890026cea9SSean Farley   MatStructure flg2;
27900026cea9SSean Farley 
27910026cea9SSean Farley   PetscFunctionBegin;
27920026cea9SSean Farley   ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
27930026cea9SSean Farley   ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr);
27940026cea9SSean Farley   ierr = MatMult(A,Xdot,F);CHKERRQ(ierr);
27950026cea9SSean Farley   PetscFunctionReturn(0);
27960026cea9SSean Farley }
27970026cea9SSean Farley 
27980026cea9SSean Farley #undef __FUNCT__
27990026cea9SSean Farley #define __FUNCT__ "TSComputeIJacobianConstant"
28000026cea9SSean Farley /*@C
280124e94211SJed Brown    TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent.
28020026cea9SSean Farley 
28030026cea9SSean Farley    Collective on TS
28040026cea9SSean Farley 
28050026cea9SSean Farley    Input Arguments:
28060026cea9SSean Farley +  ts - time stepping context
28070026cea9SSean Farley .  t - time at which to evaluate
28080026cea9SSean Farley .  X - state at which to evaluate
28090026cea9SSean Farley .  Xdot - time derivative of state vector
28100026cea9SSean Farley .  shift - shift to apply
28110026cea9SSean Farley -  ctx - context
28120026cea9SSean Farley 
28130026cea9SSean Farley    Output Arguments:
28140026cea9SSean Farley +  A - pointer to operator
28150026cea9SSean Farley .  B - pointer to preconditioning matrix
28160026cea9SSean Farley -  flg - matrix structure flag
28170026cea9SSean Farley 
28180026cea9SSean Farley    Level: intermediate
28190026cea9SSean Farley 
28200026cea9SSean Farley    Notes:
28210026cea9SSean Farley    This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems.
28220026cea9SSean Farley 
28230026cea9SSean Farley .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear()
28240026cea9SSean Farley @*/
28250026cea9SSean Farley PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx)
28260026cea9SSean Farley {
28270026cea9SSean Farley 
28280026cea9SSean Farley   PetscFunctionBegin;
28290026cea9SSean Farley   *flg = SAME_PRECONDITIONER;
28300026cea9SSean Farley   PetscFunctionReturn(0);
28310026cea9SSean Farley }
28320026cea9SSean Farley 
28334af1b03aSJed Brown 
28344af1b03aSJed Brown #undef __FUNCT__
28354af1b03aSJed Brown #define __FUNCT__ "TSGetConvergedReason"
28364af1b03aSJed Brown /*@
28374af1b03aSJed Brown    TSGetConvergedReason - Gets the reason the TS iteration was stopped.
28384af1b03aSJed Brown 
28394af1b03aSJed Brown    Not Collective
28404af1b03aSJed Brown 
28414af1b03aSJed Brown    Input Parameter:
28424af1b03aSJed Brown .  ts - the TS context
28434af1b03aSJed Brown 
28444af1b03aSJed Brown    Output Parameter:
28454af1b03aSJed Brown .  reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the
28464af1b03aSJed Brown             manual pages for the individual convergence tests for complete lists
28474af1b03aSJed Brown 
28484af1b03aSJed Brown    Level: intermediate
28494af1b03aSJed Brown 
2850cd652676SJed Brown    Notes:
2851cd652676SJed Brown    Can only be called after the call to TSSolve() is complete.
28524af1b03aSJed Brown 
28534af1b03aSJed Brown .keywords: TS, nonlinear, set, convergence, test
28544af1b03aSJed Brown 
28554af1b03aSJed Brown .seealso: TSSetConvergenceTest(), TSConvergedReason
28564af1b03aSJed Brown @*/
28574af1b03aSJed Brown PetscErrorCode  TSGetConvergedReason(TS ts,TSConvergedReason *reason)
28584af1b03aSJed Brown {
28594af1b03aSJed Brown   PetscFunctionBegin;
28604af1b03aSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
28614af1b03aSJed Brown   PetscValidPointer(reason,2);
28624af1b03aSJed Brown   *reason = ts->reason;
28634af1b03aSJed Brown   PetscFunctionReturn(0);
28644af1b03aSJed Brown }
28654af1b03aSJed Brown 
2866fb1732b5SBarry Smith #undef __FUNCT__
28675ef26d82SJed Brown #define __FUNCT__ "TSGetSNESIterations"
28689f67acb7SJed Brown /*@
28695ef26d82SJed Brown    TSGetSNESIterations - Gets the total number of nonlinear iterations
28709f67acb7SJed Brown    used by the time integrator.
28719f67acb7SJed Brown 
28729f67acb7SJed Brown    Not Collective
28739f67acb7SJed Brown 
28749f67acb7SJed Brown    Input Parameter:
28759f67acb7SJed Brown .  ts - TS context
28769f67acb7SJed Brown 
28779f67acb7SJed Brown    Output Parameter:
28789f67acb7SJed Brown .  nits - number of nonlinear iterations
28799f67acb7SJed Brown 
28809f67acb7SJed Brown    Notes:
28819f67acb7SJed Brown    This counter is reset to zero for each successive call to TSSolve().
28829f67acb7SJed Brown 
28839f67acb7SJed Brown    Level: intermediate
28849f67acb7SJed Brown 
28859f67acb7SJed Brown .keywords: TS, get, number, nonlinear, iterations
28869f67acb7SJed Brown 
28875ef26d82SJed Brown .seealso:  TSGetKSPIterations()
28889f67acb7SJed Brown @*/
28895ef26d82SJed Brown PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits)
28909f67acb7SJed Brown {
28919f67acb7SJed Brown   PetscFunctionBegin;
28929f67acb7SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
28939f67acb7SJed Brown   PetscValidIntPointer(nits,2);
28945ef26d82SJed Brown   *nits = ts->snes_its;
28959f67acb7SJed Brown   PetscFunctionReturn(0);
28969f67acb7SJed Brown }
28979f67acb7SJed Brown 
28989f67acb7SJed Brown #undef __FUNCT__
28995ef26d82SJed Brown #define __FUNCT__ "TSGetKSPIterations"
29009f67acb7SJed Brown /*@
29015ef26d82SJed Brown    TSGetKSPIterations - Gets the total number of linear iterations
29029f67acb7SJed Brown    used by the time integrator.
29039f67acb7SJed Brown 
29049f67acb7SJed Brown    Not Collective
29059f67acb7SJed Brown 
29069f67acb7SJed Brown    Input Parameter:
29079f67acb7SJed Brown .  ts - TS context
29089f67acb7SJed Brown 
29099f67acb7SJed Brown    Output Parameter:
29109f67acb7SJed Brown .  lits - number of linear iterations
29119f67acb7SJed Brown 
29129f67acb7SJed Brown    Notes:
29139f67acb7SJed Brown    This counter is reset to zero for each successive call to TSSolve().
29149f67acb7SJed Brown 
29159f67acb7SJed Brown    Level: intermediate
29169f67acb7SJed Brown 
29179f67acb7SJed Brown .keywords: TS, get, number, linear, iterations
29189f67acb7SJed Brown 
29195ef26d82SJed Brown .seealso:  TSGetSNESIterations(), SNESGetKSPIterations()
29209f67acb7SJed Brown @*/
29215ef26d82SJed Brown PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits)
29229f67acb7SJed Brown {
29239f67acb7SJed Brown   PetscFunctionBegin;
29249f67acb7SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
29259f67acb7SJed Brown   PetscValidIntPointer(lits,2);
29265ef26d82SJed Brown   *lits = ts->ksp_its;
29279f67acb7SJed Brown   PetscFunctionReturn(0);
29289f67acb7SJed Brown }
29299f67acb7SJed Brown 
29309f67acb7SJed Brown #undef __FUNCT__
2931cef5090cSJed Brown #define __FUNCT__ "TSGetStepRejections"
2932cef5090cSJed Brown /*@
2933cef5090cSJed Brown    TSGetStepRejections - Gets the total number of rejected steps.
2934cef5090cSJed Brown 
2935cef5090cSJed Brown    Not Collective
2936cef5090cSJed Brown 
2937cef5090cSJed Brown    Input Parameter:
2938cef5090cSJed Brown .  ts - TS context
2939cef5090cSJed Brown 
2940cef5090cSJed Brown    Output Parameter:
2941cef5090cSJed Brown .  rejects - number of steps rejected
2942cef5090cSJed Brown 
2943cef5090cSJed Brown    Notes:
2944cef5090cSJed Brown    This counter is reset to zero for each successive call to TSSolve().
2945cef5090cSJed Brown 
2946cef5090cSJed Brown    Level: intermediate
2947cef5090cSJed Brown 
2948cef5090cSJed Brown .keywords: TS, get, number
2949cef5090cSJed Brown 
29505ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails()
2951cef5090cSJed Brown @*/
2952cef5090cSJed Brown PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects)
2953cef5090cSJed Brown {
2954cef5090cSJed Brown   PetscFunctionBegin;
2955cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2956cef5090cSJed Brown   PetscValidIntPointer(rejects,2);
2957cef5090cSJed Brown   *rejects = ts->reject;
2958cef5090cSJed Brown   PetscFunctionReturn(0);
2959cef5090cSJed Brown }
2960cef5090cSJed Brown 
2961cef5090cSJed Brown #undef __FUNCT__
2962cef5090cSJed Brown #define __FUNCT__ "TSGetSNESFailures"
2963cef5090cSJed Brown /*@
2964cef5090cSJed Brown    TSGetSNESFailures - Gets the total number of failed SNES solves
2965cef5090cSJed Brown 
2966cef5090cSJed Brown    Not Collective
2967cef5090cSJed Brown 
2968cef5090cSJed Brown    Input Parameter:
2969cef5090cSJed Brown .  ts - TS context
2970cef5090cSJed Brown 
2971cef5090cSJed Brown    Output Parameter:
2972cef5090cSJed Brown .  fails - number of failed nonlinear solves
2973cef5090cSJed Brown 
2974cef5090cSJed Brown    Notes:
2975cef5090cSJed Brown    This counter is reset to zero for each successive call to TSSolve().
2976cef5090cSJed Brown 
2977cef5090cSJed Brown    Level: intermediate
2978cef5090cSJed Brown 
2979cef5090cSJed Brown .keywords: TS, get, number
2980cef5090cSJed Brown 
29815ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures()
2982cef5090cSJed Brown @*/
2983cef5090cSJed Brown PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails)
2984cef5090cSJed Brown {
2985cef5090cSJed Brown   PetscFunctionBegin;
2986cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2987cef5090cSJed Brown   PetscValidIntPointer(fails,2);
2988cef5090cSJed Brown   *fails = ts->num_snes_failures;
2989cef5090cSJed Brown   PetscFunctionReturn(0);
2990cef5090cSJed Brown }
2991cef5090cSJed Brown 
2992cef5090cSJed Brown #undef __FUNCT__
2993cef5090cSJed Brown #define __FUNCT__ "TSSetMaxStepRejections"
2994cef5090cSJed Brown /*@
2995cef5090cSJed Brown    TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails
2996cef5090cSJed Brown 
2997cef5090cSJed Brown    Not Collective
2998cef5090cSJed Brown 
2999cef5090cSJed Brown    Input Parameter:
3000cef5090cSJed Brown +  ts - TS context
3001cef5090cSJed Brown -  rejects - maximum number of rejected steps, pass -1 for unlimited
3002cef5090cSJed Brown 
3003cef5090cSJed Brown    Notes:
3004cef5090cSJed Brown    The counter is reset to zero for each step
3005cef5090cSJed Brown 
3006cef5090cSJed Brown    Options Database Key:
3007cef5090cSJed Brown  .  -ts_max_reject - Maximum number of step rejections before a step fails
3008cef5090cSJed Brown 
3009cef5090cSJed Brown    Level: intermediate
3010cef5090cSJed Brown 
3011cef5090cSJed Brown .keywords: TS, set, maximum, number
3012cef5090cSJed Brown 
30135ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3014cef5090cSJed Brown @*/
3015cef5090cSJed Brown PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects)
3016cef5090cSJed Brown {
3017cef5090cSJed Brown   PetscFunctionBegin;
3018cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3019cef5090cSJed Brown   ts->max_reject = rejects;
3020cef5090cSJed Brown   PetscFunctionReturn(0);
3021cef5090cSJed Brown }
3022cef5090cSJed Brown 
3023cef5090cSJed Brown #undef __FUNCT__
3024cef5090cSJed Brown #define __FUNCT__ "TSSetMaxSNESFailures"
3025cef5090cSJed Brown /*@
3026cef5090cSJed Brown    TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves
3027cef5090cSJed Brown 
3028cef5090cSJed Brown    Not Collective
3029cef5090cSJed Brown 
3030cef5090cSJed Brown    Input Parameter:
3031cef5090cSJed Brown +  ts - TS context
3032cef5090cSJed Brown -  fails - maximum number of failed nonlinear solves, pass -1 for unlimited
3033cef5090cSJed Brown 
3034cef5090cSJed Brown    Notes:
3035cef5090cSJed Brown    The counter is reset to zero for each successive call to TSSolve().
3036cef5090cSJed Brown 
3037cef5090cSJed Brown    Options Database Key:
3038cef5090cSJed Brown  .  -ts_max_snes_failures - Maximum number of nonlinear solve failures
3039cef5090cSJed Brown 
3040cef5090cSJed Brown    Level: intermediate
3041cef5090cSJed Brown 
3042cef5090cSJed Brown .keywords: TS, set, maximum, number
3043cef5090cSJed Brown 
30445ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason()
3045cef5090cSJed Brown @*/
3046cef5090cSJed Brown PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails)
3047cef5090cSJed Brown {
3048cef5090cSJed Brown   PetscFunctionBegin;
3049cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3050cef5090cSJed Brown   ts->max_snes_failures = fails;
3051cef5090cSJed Brown   PetscFunctionReturn(0);
3052cef5090cSJed Brown }
3053cef5090cSJed Brown 
3054cef5090cSJed Brown #undef __FUNCT__
3055cef5090cSJed Brown #define __FUNCT__ "TSSetErrorIfStepFails()"
3056cef5090cSJed Brown /*@
3057cef5090cSJed Brown    TSSetErrorIfStepFails - Error if no step succeeds
3058cef5090cSJed Brown 
3059cef5090cSJed Brown    Not Collective
3060cef5090cSJed Brown 
3061cef5090cSJed Brown    Input Parameter:
3062cef5090cSJed Brown +  ts - TS context
3063cef5090cSJed Brown -  err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure
3064cef5090cSJed Brown 
3065cef5090cSJed Brown    Options Database Key:
3066cef5090cSJed Brown  .  -ts_error_if_step_fails - Error if no step succeeds
3067cef5090cSJed Brown 
3068cef5090cSJed Brown    Level: intermediate
3069cef5090cSJed Brown 
3070cef5090cSJed Brown .keywords: TS, set, error
3071cef5090cSJed Brown 
30725ef26d82SJed Brown .seealso:  TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason()
3073cef5090cSJed Brown @*/
3074cef5090cSJed Brown PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err)
3075cef5090cSJed Brown {
3076cef5090cSJed Brown   PetscFunctionBegin;
3077cef5090cSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3078cef5090cSJed Brown   ts->errorifstepfailed = err;
3079cef5090cSJed Brown   PetscFunctionReturn(0);
3080cef5090cSJed Brown }
3081cef5090cSJed Brown 
3082cef5090cSJed Brown #undef __FUNCT__
3083fb1732b5SBarry Smith #define __FUNCT__ "TSMonitorSolutionBinary"
3084fb1732b5SBarry Smith /*@C
3085fb1732b5SBarry Smith    TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file
3086fb1732b5SBarry Smith 
3087fb1732b5SBarry Smith    Collective on TS
3088fb1732b5SBarry Smith 
3089fb1732b5SBarry Smith    Input Parameters:
3090fb1732b5SBarry Smith +  ts - the TS context
3091fb1732b5SBarry Smith .  step - current time-step
3092fb1732b5SBarry Smith .  ptime - current time
3093ed81e22dSJed Brown .  x - current state
3094fb1732b5SBarry Smith -  viewer - binary viewer
3095fb1732b5SBarry Smith 
3096fb1732b5SBarry Smith    Level: intermediate
3097fb1732b5SBarry Smith 
3098fb1732b5SBarry Smith .keywords: TS,  vector, monitor, view
3099fb1732b5SBarry Smith 
3100fb1732b5SBarry Smith .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3101fb1732b5SBarry Smith @*/
3102ed81e22dSJed Brown PetscErrorCode  TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer)
3103fb1732b5SBarry Smith {
3104fb1732b5SBarry Smith   PetscErrorCode       ierr;
3105ed81e22dSJed Brown   PetscViewer          v = (PetscViewer)viewer;
3106fb1732b5SBarry Smith 
3107fb1732b5SBarry Smith   PetscFunctionBegin;
3108ed81e22dSJed Brown   ierr = VecView(x,v);CHKERRQ(ierr);
3109ed81e22dSJed Brown   PetscFunctionReturn(0);
3110ed81e22dSJed Brown }
3111ed81e22dSJed Brown 
3112ed81e22dSJed Brown #undef __FUNCT__
3113ed81e22dSJed Brown #define __FUNCT__ "TSMonitorSolutionVTK"
3114ed81e22dSJed Brown /*@C
3115ed81e22dSJed Brown    TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep.
3116ed81e22dSJed Brown 
3117ed81e22dSJed Brown    Collective on TS
3118ed81e22dSJed Brown 
3119ed81e22dSJed Brown    Input Parameters:
3120ed81e22dSJed Brown +  ts - the TS context
3121ed81e22dSJed Brown .  step - current time-step
3122ed81e22dSJed Brown .  ptime - current time
3123ed81e22dSJed Brown .  x - current state
3124ed81e22dSJed Brown -  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3125ed81e22dSJed Brown 
3126ed81e22dSJed Brown    Level: intermediate
3127ed81e22dSJed Brown 
3128ed81e22dSJed Brown    Notes:
3129ed81e22dSJed 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.
3130ed81e22dSJed Brown    These are named according to the file name template.
3131ed81e22dSJed Brown 
3132ed81e22dSJed Brown    This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy().
3133ed81e22dSJed Brown 
3134ed81e22dSJed Brown .keywords: TS,  vector, monitor, view
3135ed81e22dSJed Brown 
3136ed81e22dSJed Brown .seealso: TSMonitorSet(), TSMonitorDefault(), VecView()
3137ed81e22dSJed Brown @*/
3138ed81e22dSJed Brown PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate)
3139ed81e22dSJed Brown {
3140ed81e22dSJed Brown   PetscErrorCode ierr;
3141ed81e22dSJed Brown   char           filename[PETSC_MAX_PATH_LEN];
3142ed81e22dSJed Brown   PetscViewer    viewer;
3143ed81e22dSJed Brown 
3144ed81e22dSJed Brown   PetscFunctionBegin;
3145ed81e22dSJed Brown   ierr = PetscSNPrintf(filename,sizeof filename,(const char*)filenametemplate,step);CHKERRQ(ierr);
3146ed81e22dSJed Brown   ierr = PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
3147fb1732b5SBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
3148ed81e22dSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3149ed81e22dSJed Brown   PetscFunctionReturn(0);
3150ed81e22dSJed Brown }
3151ed81e22dSJed Brown 
3152ed81e22dSJed Brown #undef __FUNCT__
3153ed81e22dSJed Brown #define __FUNCT__ "TSMonitorSolutionVTKDestroy"
3154ed81e22dSJed Brown /*@C
3155ed81e22dSJed Brown    TSMonitorSolutionVTKDestroy - Destroy context for monitoring
3156ed81e22dSJed Brown 
3157ed81e22dSJed Brown    Collective on TS
3158ed81e22dSJed Brown 
3159ed81e22dSJed Brown    Input Parameters:
3160ed81e22dSJed Brown .  filenametemplate - string containing a format specifier for the integer time step (e.g. %03D)
3161ed81e22dSJed Brown 
3162ed81e22dSJed Brown    Level: intermediate
3163ed81e22dSJed Brown 
3164ed81e22dSJed Brown    Note:
3165ed81e22dSJed Brown    This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK().
3166ed81e22dSJed Brown 
3167ed81e22dSJed Brown .keywords: TS,  vector, monitor, view
3168ed81e22dSJed Brown 
3169ed81e22dSJed Brown .seealso: TSMonitorSet(), TSMonitorSolutionVTK()
3170ed81e22dSJed Brown @*/
3171ed81e22dSJed Brown PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate)
3172ed81e22dSJed Brown {
3173ed81e22dSJed Brown   PetscErrorCode ierr;
3174ed81e22dSJed Brown 
3175ed81e22dSJed Brown   PetscFunctionBegin;
3176ed81e22dSJed Brown   ierr = PetscFree(*(char**)filenametemplate);CHKERRQ(ierr);
3177fb1732b5SBarry Smith   PetscFunctionReturn(0);
3178fb1732b5SBarry Smith }
3179fb1732b5SBarry Smith 
318084df9cb4SJed Brown #undef __FUNCT__
318184df9cb4SJed Brown #define __FUNCT__ "TSGetAdapt"
318284df9cb4SJed Brown /*@
318384df9cb4SJed Brown    TSGetAdapt - Get the adaptive controller context for the current method
318484df9cb4SJed Brown 
3185ed81e22dSJed Brown    Collective on TS if controller has not been created yet
318684df9cb4SJed Brown 
318784df9cb4SJed Brown    Input Arguments:
3188ed81e22dSJed Brown .  ts - time stepping context
318984df9cb4SJed Brown 
319084df9cb4SJed Brown    Output Arguments:
3191ed81e22dSJed Brown .  adapt - adaptive controller
319284df9cb4SJed Brown 
319384df9cb4SJed Brown    Level: intermediate
319484df9cb4SJed Brown 
3195ed81e22dSJed Brown .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose()
319684df9cb4SJed Brown @*/
319784df9cb4SJed Brown PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt)
319884df9cb4SJed Brown {
319984df9cb4SJed Brown   PetscErrorCode ierr;
320084df9cb4SJed Brown 
320184df9cb4SJed Brown   PetscFunctionBegin;
320284df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
320384df9cb4SJed Brown   PetscValidPointer(adapt,2);
320484df9cb4SJed Brown   if (!ts->adapt) {
320584df9cb4SJed Brown     ierr = TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);CHKERRQ(ierr);
32061c3436cfSJed Brown     ierr = PetscLogObjectParent(ts,ts->adapt);CHKERRQ(ierr);
32071c3436cfSJed Brown     ierr = PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);CHKERRQ(ierr);
320884df9cb4SJed Brown   }
320984df9cb4SJed Brown   *adapt = ts->adapt;
321084df9cb4SJed Brown   PetscFunctionReturn(0);
321184df9cb4SJed Brown }
3212d6ebe24aSShri Abhyankar 
3213d6ebe24aSShri Abhyankar #undef __FUNCT__
32141c3436cfSJed Brown #define __FUNCT__ "TSSetTolerances"
32151c3436cfSJed Brown /*@
32161c3436cfSJed Brown    TSSetTolerances - Set tolerances for local truncation error when using adaptive controller
32171c3436cfSJed Brown 
32181c3436cfSJed Brown    Logically Collective
32191c3436cfSJed Brown 
32201c3436cfSJed Brown    Input Arguments:
32211c3436cfSJed Brown +  ts - time integration context
32221c3436cfSJed Brown .  atol - scalar absolute tolerances, PETSC_DECIDE to leave current value
32231c3436cfSJed Brown .  vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present
32241c3436cfSJed Brown .  rtol - scalar relative tolerances, PETSC_DECIDE to leave current value
32251c3436cfSJed Brown -  vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present
32261c3436cfSJed Brown 
32271c3436cfSJed Brown    Level: beginner
32281c3436cfSJed Brown 
3229c5033834SJed Brown .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances()
32301c3436cfSJed Brown @*/
32311c3436cfSJed Brown PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol)
32321c3436cfSJed Brown {
32331c3436cfSJed Brown   PetscErrorCode ierr;
32341c3436cfSJed Brown 
32351c3436cfSJed Brown   PetscFunctionBegin;
3236c5033834SJed Brown   if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol;
32371c3436cfSJed Brown   if (vatol) {
32381c3436cfSJed Brown     ierr = PetscObjectReference((PetscObject)vatol);CHKERRQ(ierr);
32391c3436cfSJed Brown     ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr);
32401c3436cfSJed Brown     ts->vatol = vatol;
32411c3436cfSJed Brown   }
3242c5033834SJed Brown   if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol;
32431c3436cfSJed Brown   if (vrtol) {
32441c3436cfSJed Brown     ierr = PetscObjectReference((PetscObject)vrtol);CHKERRQ(ierr);
32451c3436cfSJed Brown     ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr);
32461c3436cfSJed Brown     ts->vrtol = vrtol;
32471c3436cfSJed Brown   }
32481c3436cfSJed Brown   PetscFunctionReturn(0);
32491c3436cfSJed Brown }
32501c3436cfSJed Brown 
32511c3436cfSJed Brown #undef __FUNCT__
3252c5033834SJed Brown #define __FUNCT__ "TSGetTolerances"
3253c5033834SJed Brown /*@
3254c5033834SJed Brown    TSGetTolerances - Get tolerances for local truncation error when using adaptive controller
3255c5033834SJed Brown 
3256c5033834SJed Brown    Logically Collective
3257c5033834SJed Brown 
3258c5033834SJed Brown    Input Arguments:
3259c5033834SJed Brown .  ts - time integration context
3260c5033834SJed Brown 
3261c5033834SJed Brown    Output Arguments:
3262c5033834SJed Brown +  atol - scalar absolute tolerances, PETSC_NULL to ignore
3263c5033834SJed Brown .  vatol - vector of absolute tolerances, PETSC_NULL to ignore
3264c5033834SJed Brown .  rtol - scalar relative tolerances, PETSC_NULL to ignore
3265c5033834SJed Brown -  vrtol - vector of relative tolerances, PETSC_NULL to ignore
3266c5033834SJed Brown 
3267c5033834SJed Brown    Level: beginner
3268c5033834SJed Brown 
3269c5033834SJed Brown .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances()
3270c5033834SJed Brown @*/
3271c5033834SJed Brown PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol)
3272c5033834SJed Brown {
3273c5033834SJed Brown 
3274c5033834SJed Brown   PetscFunctionBegin;
3275c5033834SJed Brown   if (atol)  *atol  = ts->atol;
3276c5033834SJed Brown   if (vatol) *vatol = ts->vatol;
3277c5033834SJed Brown   if (rtol)  *rtol  = ts->rtol;
3278c5033834SJed Brown   if (vrtol) *vrtol = ts->vrtol;
3279c5033834SJed Brown   PetscFunctionReturn(0);
3280c5033834SJed Brown }
3281c5033834SJed Brown 
3282c5033834SJed Brown #undef __FUNCT__
32831c3436cfSJed Brown #define __FUNCT__ "TSErrorNormWRMS"
32841c3436cfSJed Brown /*@
32851c3436cfSJed Brown    TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state
32861c3436cfSJed Brown 
32871c3436cfSJed Brown    Collective on TS
32881c3436cfSJed Brown 
32891c3436cfSJed Brown    Input Arguments:
32901c3436cfSJed Brown +  ts - time stepping context
32911c3436cfSJed Brown -  Y - state vector to be compared to ts->vec_sol
32921c3436cfSJed Brown 
32931c3436cfSJed Brown    Output Arguments:
32941c3436cfSJed Brown .  norm - weighted norm, a value of 1.0 is considered small
32951c3436cfSJed Brown 
32961c3436cfSJed Brown    Level: developer
32971c3436cfSJed Brown 
32981c3436cfSJed Brown .seealso: TSSetTolerances()
32991c3436cfSJed Brown @*/
33001c3436cfSJed Brown PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm)
33011c3436cfSJed Brown {
33021c3436cfSJed Brown   PetscErrorCode ierr;
33031c3436cfSJed Brown   PetscInt i,n,N;
33041c3436cfSJed Brown   const PetscScalar *x,*y;
33051c3436cfSJed Brown   Vec X;
33061c3436cfSJed Brown   PetscReal sum,gsum;
33071c3436cfSJed Brown 
33081c3436cfSJed Brown   PetscFunctionBegin;
33091c3436cfSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
33101c3436cfSJed Brown   PetscValidHeaderSpecific(Y,VEC_CLASSID,2);
33111c3436cfSJed Brown   PetscValidPointer(norm,3);
33121c3436cfSJed Brown   X = ts->vec_sol;
33131c3436cfSJed Brown   PetscCheckSameTypeAndComm(X,1,Y,2);
33141c3436cfSJed Brown   if (X == Y) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector");
33151c3436cfSJed Brown 
33161c3436cfSJed Brown   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
33171c3436cfSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
33181c3436cfSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
33191c3436cfSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
33201c3436cfSJed Brown   sum = 0.;
3321ceefc113SJed Brown   if (ts->vatol && ts->vrtol) {
3322ceefc113SJed Brown     const PetscScalar *atol,*rtol;
3323ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3324ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3325ceefc113SJed Brown     for (i=0; i<n; i++) {
3326ceefc113SJed Brown       PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3327ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3328ceefc113SJed Brown     }
3329ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3330ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3331ceefc113SJed Brown   } else if (ts->vatol) {       /* vector atol, scalar rtol */
3332ceefc113SJed Brown     const PetscScalar *atol;
3333ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3334ceefc113SJed Brown     for (i=0; i<n; i++) {
3335ceefc113SJed Brown       PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3336ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3337ceefc113SJed Brown     }
3338ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr);
3339ceefc113SJed Brown   } else if (ts->vrtol) {       /* scalar atol, vector rtol */
3340ceefc113SJed Brown     const PetscScalar *rtol;
3341ceefc113SJed Brown     ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3342ceefc113SJed Brown     for (i=0; i<n; i++) {
3343ceefc113SJed Brown       PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
3344ceefc113SJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
3345ceefc113SJed Brown     }
3346ceefc113SJed Brown     ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr);
3347ceefc113SJed Brown   } else {                      /* scalar atol, scalar rtol */
33481c3436cfSJed Brown     for (i=0; i<n; i++) {
33491c3436cfSJed Brown       PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i]));
33501c3436cfSJed Brown       sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol);
33511c3436cfSJed Brown     }
3352ceefc113SJed Brown   }
33531c3436cfSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
33541c3436cfSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
33551c3436cfSJed Brown 
33561c3436cfSJed Brown   ierr = MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);CHKERRQ(ierr);
33571c3436cfSJed Brown   *norm = PetscSqrtReal(gsum / N);
33581c3436cfSJed Brown   if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
33591c3436cfSJed Brown   PetscFunctionReturn(0);
33601c3436cfSJed Brown }
33611c3436cfSJed Brown 
33621c3436cfSJed Brown #undef __FUNCT__
33638d59e960SJed Brown #define __FUNCT__ "TSSetCFLTimeLocal"
33648d59e960SJed Brown /*@
33658d59e960SJed Brown    TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler
33668d59e960SJed Brown 
33678d59e960SJed Brown    Logically Collective on TS
33688d59e960SJed Brown 
33698d59e960SJed Brown    Input Arguments:
33708d59e960SJed Brown +  ts - time stepping context
33718d59e960SJed Brown -  cfltime - maximum stable time step if using forward Euler (value can be different on each process)
33728d59e960SJed Brown 
33738d59e960SJed Brown    Note:
33748d59e960SJed Brown    After calling this function, the global CFL time can be obtained by calling TSGetCFLTime()
33758d59e960SJed Brown 
33768d59e960SJed Brown    Level: intermediate
33778d59e960SJed Brown 
33788d59e960SJed Brown .seealso: TSGetCFLTime(), TSADAPTCFL
33798d59e960SJed Brown @*/
33808d59e960SJed Brown PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime)
33818d59e960SJed Brown {
33828d59e960SJed Brown 
33838d59e960SJed Brown   PetscFunctionBegin;
33848d59e960SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
33858d59e960SJed Brown   ts->cfltime_local = cfltime;
33868d59e960SJed Brown   ts->cfltime = -1.;
33878d59e960SJed Brown   PetscFunctionReturn(0);
33888d59e960SJed Brown }
33898d59e960SJed Brown 
33908d59e960SJed Brown #undef __FUNCT__
33918d59e960SJed Brown #define __FUNCT__ "TSGetCFLTime"
33928d59e960SJed Brown /*@
33938d59e960SJed Brown    TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler
33948d59e960SJed Brown 
33958d59e960SJed Brown    Collective on TS
33968d59e960SJed Brown 
33978d59e960SJed Brown    Input Arguments:
33988d59e960SJed Brown .  ts - time stepping context
33998d59e960SJed Brown 
34008d59e960SJed Brown    Output Arguments:
34018d59e960SJed Brown .  cfltime - maximum stable time step for forward Euler
34028d59e960SJed Brown 
34038d59e960SJed Brown    Level: advanced
34048d59e960SJed Brown 
34058d59e960SJed Brown .seealso: TSSetCFLTimeLocal()
34068d59e960SJed Brown @*/
34078d59e960SJed Brown PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime)
34088d59e960SJed Brown {
34098d59e960SJed Brown   PetscErrorCode ierr;
34108d59e960SJed Brown 
34118d59e960SJed Brown   PetscFunctionBegin;
34128d59e960SJed Brown   if (ts->cfltime < 0) {
34138d59e960SJed Brown     ierr = MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr);
34148d59e960SJed Brown   }
34158d59e960SJed Brown   *cfltime = ts->cfltime;
34168d59e960SJed Brown   PetscFunctionReturn(0);
34178d59e960SJed Brown }
34188d59e960SJed Brown 
34198d59e960SJed Brown #undef __FUNCT__
3420d6ebe24aSShri Abhyankar #define __FUNCT__ "TSVISetVariableBounds"
3421d6ebe24aSShri Abhyankar /*@
3422d6ebe24aSShri Abhyankar    TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu
3423d6ebe24aSShri Abhyankar 
3424d6ebe24aSShri Abhyankar    Input Parameters:
3425d6ebe24aSShri Abhyankar .  ts   - the TS context.
3426d6ebe24aSShri Abhyankar .  xl   - lower bound.
3427d6ebe24aSShri Abhyankar .  xu   - upper bound.
3428d6ebe24aSShri Abhyankar 
3429d6ebe24aSShri Abhyankar    Notes:
3430d6ebe24aSShri Abhyankar    If this routine is not called then the lower and upper bounds are set to
343133c0c2ddSJed Brown    SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp().
3432d6ebe24aSShri Abhyankar 
34332bd2b0e6SSatish Balay    Level: advanced
34342bd2b0e6SSatish Balay 
3435d6ebe24aSShri Abhyankar @*/
3436d6ebe24aSShri Abhyankar PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu)
3437d6ebe24aSShri Abhyankar {
3438d6ebe24aSShri Abhyankar   PetscErrorCode ierr;
3439d6ebe24aSShri Abhyankar   SNES           snes;
3440d6ebe24aSShri Abhyankar 
3441d6ebe24aSShri Abhyankar   PetscFunctionBegin;
3442d6ebe24aSShri Abhyankar   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3443d6ebe24aSShri Abhyankar   ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr);
3444d6ebe24aSShri Abhyankar   PetscFunctionReturn(0);
3445d6ebe24aSShri Abhyankar }
3446d6ebe24aSShri Abhyankar 
3447325fc9f4SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
3448c6db04a5SJed Brown #include <mex.h>
3449325fc9f4SBarry Smith 
3450325fc9f4SBarry Smith typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext;
3451325fc9f4SBarry Smith 
3452325fc9f4SBarry Smith #undef __FUNCT__
3453325fc9f4SBarry Smith #define __FUNCT__ "TSComputeFunction_Matlab"
3454325fc9f4SBarry Smith /*
3455325fc9f4SBarry Smith    TSComputeFunction_Matlab - Calls the function that has been set with
3456325fc9f4SBarry Smith                          TSSetFunctionMatlab().
3457325fc9f4SBarry Smith 
3458325fc9f4SBarry Smith    Collective on TS
3459325fc9f4SBarry Smith 
3460325fc9f4SBarry Smith    Input Parameters:
3461325fc9f4SBarry Smith +  snes - the TS context
3462325fc9f4SBarry Smith -  x - input vector
3463325fc9f4SBarry Smith 
3464325fc9f4SBarry Smith    Output Parameter:
3465325fc9f4SBarry Smith .  y - function vector, as set by TSSetFunction()
3466325fc9f4SBarry Smith 
3467325fc9f4SBarry Smith    Notes:
3468325fc9f4SBarry Smith    TSComputeFunction() is typically used within nonlinear solvers
3469325fc9f4SBarry Smith    implementations, so most users would not generally call this routine
3470325fc9f4SBarry Smith    themselves.
3471325fc9f4SBarry Smith 
3472325fc9f4SBarry Smith    Level: developer
3473325fc9f4SBarry Smith 
3474325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
3475325fc9f4SBarry Smith 
3476325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3477325fc9f4SBarry Smith */
34787087cfbeSBarry Smith PetscErrorCode  TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx)
3479325fc9f4SBarry Smith {
3480325fc9f4SBarry Smith   PetscErrorCode   ierr;
3481325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3482325fc9f4SBarry Smith   int              nlhs = 1,nrhs = 7;
3483325fc9f4SBarry Smith   mxArray          *plhs[1],*prhs[7];
3484325fc9f4SBarry Smith   long long int    lx = 0,lxdot = 0,ly = 0,ls = 0;
3485325fc9f4SBarry Smith 
3486325fc9f4SBarry Smith   PetscFunctionBegin;
3487325fc9f4SBarry Smith   PetscValidHeaderSpecific(snes,TS_CLASSID,1);
3488325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3489325fc9f4SBarry Smith   PetscValidHeaderSpecific(xdot,VEC_CLASSID,4);
3490325fc9f4SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,5);
3491325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,x,3);
3492325fc9f4SBarry Smith   PetscCheckSameComm(snes,1,y,5);
3493325fc9f4SBarry Smith 
3494325fc9f4SBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3495325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3496880f3077SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr);
3497325fc9f4SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
3498325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3499325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar(time);
3500325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
3501325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3502325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)ly);
3503325fc9f4SBarry Smith   prhs[5] =  mxCreateString(sctx->funcname);
3504325fc9f4SBarry Smith   prhs[6] =  sctx->ctx;
3505325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr);
3506325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3507325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
3508325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
3509325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
3510325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
3511325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
3512325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
3513325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
3514325fc9f4SBarry Smith   PetscFunctionReturn(0);
3515325fc9f4SBarry Smith }
3516325fc9f4SBarry Smith 
3517325fc9f4SBarry Smith 
3518325fc9f4SBarry Smith #undef __FUNCT__
3519325fc9f4SBarry Smith #define __FUNCT__ "TSSetFunctionMatlab"
3520325fc9f4SBarry Smith /*
3521325fc9f4SBarry Smith    TSSetFunctionMatlab - Sets the function evaluation routine and function
3522325fc9f4SBarry Smith    vector for use by the TS routines in solving ODEs
3523e3c5b3baSBarry Smith    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
3524325fc9f4SBarry Smith 
3525325fc9f4SBarry Smith    Logically Collective on TS
3526325fc9f4SBarry Smith 
3527325fc9f4SBarry Smith    Input Parameters:
3528325fc9f4SBarry Smith +  ts - the TS context
3529325fc9f4SBarry Smith -  func - function evaluation routine
3530325fc9f4SBarry Smith 
3531325fc9f4SBarry Smith    Calling sequence of func:
3532325fc9f4SBarry Smith $    func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx);
3533325fc9f4SBarry Smith 
3534325fc9f4SBarry Smith    Level: beginner
3535325fc9f4SBarry Smith 
3536325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
3537325fc9f4SBarry Smith 
3538325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3539325fc9f4SBarry Smith */
3540cdcf91faSSean Farley PetscErrorCode  TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx)
3541325fc9f4SBarry Smith {
3542325fc9f4SBarry Smith   PetscErrorCode  ierr;
3543325fc9f4SBarry Smith   TSMatlabContext *sctx;
3544325fc9f4SBarry Smith 
3545325fc9f4SBarry Smith   PetscFunctionBegin;
3546325fc9f4SBarry Smith   /* currently sctx is memory bleed */
3547325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3548325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3549325fc9f4SBarry Smith   /*
3550325fc9f4SBarry Smith      This should work, but it doesn't
3551325fc9f4SBarry Smith   sctx->ctx = ctx;
3552325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
3553325fc9f4SBarry Smith   */
3554325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
3555cdcf91faSSean Farley   ierr = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr);
3556325fc9f4SBarry Smith   PetscFunctionReturn(0);
3557325fc9f4SBarry Smith }
3558325fc9f4SBarry Smith 
3559325fc9f4SBarry Smith #undef __FUNCT__
3560325fc9f4SBarry Smith #define __FUNCT__ "TSComputeJacobian_Matlab"
3561325fc9f4SBarry Smith /*
3562325fc9f4SBarry Smith    TSComputeJacobian_Matlab - Calls the function that has been set with
3563325fc9f4SBarry Smith                          TSSetJacobianMatlab().
3564325fc9f4SBarry Smith 
3565325fc9f4SBarry Smith    Collective on TS
3566325fc9f4SBarry Smith 
3567325fc9f4SBarry Smith    Input Parameters:
3568cdcf91faSSean Farley +  ts - the TS context
3569325fc9f4SBarry Smith .  x - input vector
3570325fc9f4SBarry Smith .  A, B - the matrices
3571325fc9f4SBarry Smith -  ctx - user context
3572325fc9f4SBarry Smith 
3573325fc9f4SBarry Smith    Output Parameter:
3574325fc9f4SBarry Smith .  flag - structure of the matrix
3575325fc9f4SBarry Smith 
3576325fc9f4SBarry Smith    Level: developer
3577325fc9f4SBarry Smith 
3578325fc9f4SBarry Smith .keywords: TS, nonlinear, compute, function
3579325fc9f4SBarry Smith 
3580325fc9f4SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3581325fc9f4SBarry Smith @*/
3582cdcf91faSSean Farley PetscErrorCode  TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx)
3583325fc9f4SBarry Smith {
3584325fc9f4SBarry Smith   PetscErrorCode  ierr;
3585325fc9f4SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3586325fc9f4SBarry Smith   int             nlhs = 2,nrhs = 9;
3587325fc9f4SBarry Smith   mxArray         *plhs[2],*prhs[9];
3588325fc9f4SBarry Smith   long long int   lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0;
3589325fc9f4SBarry Smith 
3590325fc9f4SBarry Smith   PetscFunctionBegin;
3591cdcf91faSSean Farley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3592325fc9f4SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3593325fc9f4SBarry Smith 
3594325fc9f4SBarry Smith   /* call Matlab function in ctx with arguments x and y */
3595325fc9f4SBarry Smith 
3596cdcf91faSSean Farley   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
3597325fc9f4SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3598325fc9f4SBarry Smith   ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr);
3599325fc9f4SBarry Smith   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
3600325fc9f4SBarry Smith   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
3601325fc9f4SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3602325fc9f4SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)time);
3603325fc9f4SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lx);
3604325fc9f4SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lxdot);
3605325fc9f4SBarry Smith   prhs[4] =  mxCreateDoubleScalar((double)shift);
3606325fc9f4SBarry Smith   prhs[5] =  mxCreateDoubleScalar((double)lA);
3607325fc9f4SBarry Smith   prhs[6] =  mxCreateDoubleScalar((double)lB);
3608325fc9f4SBarry Smith   prhs[7] =  mxCreateString(sctx->funcname);
3609325fc9f4SBarry Smith   prhs[8] =  sctx->ctx;
3610325fc9f4SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr);
3611325fc9f4SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3612325fc9f4SBarry Smith   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
3613325fc9f4SBarry Smith   mxDestroyArray(prhs[0]);
3614325fc9f4SBarry Smith   mxDestroyArray(prhs[1]);
3615325fc9f4SBarry Smith   mxDestroyArray(prhs[2]);
3616325fc9f4SBarry Smith   mxDestroyArray(prhs[3]);
3617325fc9f4SBarry Smith   mxDestroyArray(prhs[4]);
3618325fc9f4SBarry Smith   mxDestroyArray(prhs[5]);
3619325fc9f4SBarry Smith   mxDestroyArray(prhs[6]);
3620325fc9f4SBarry Smith   mxDestroyArray(prhs[7]);
3621325fc9f4SBarry Smith   mxDestroyArray(plhs[0]);
3622325fc9f4SBarry Smith   mxDestroyArray(plhs[1]);
3623325fc9f4SBarry Smith   PetscFunctionReturn(0);
3624325fc9f4SBarry Smith }
3625325fc9f4SBarry Smith 
3626325fc9f4SBarry Smith 
3627325fc9f4SBarry Smith #undef __FUNCT__
3628325fc9f4SBarry Smith #define __FUNCT__ "TSSetJacobianMatlab"
3629325fc9f4SBarry Smith /*
3630325fc9f4SBarry Smith    TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
3631e3c5b3baSBarry 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
3632325fc9f4SBarry Smith 
3633325fc9f4SBarry Smith    Logically Collective on TS
3634325fc9f4SBarry Smith 
3635325fc9f4SBarry Smith    Input Parameters:
3636cdcf91faSSean Farley +  ts - the TS context
3637325fc9f4SBarry Smith .  A,B - Jacobian matrices
3638325fc9f4SBarry Smith .  func - function evaluation routine
3639325fc9f4SBarry Smith -  ctx - user context
3640325fc9f4SBarry Smith 
3641325fc9f4SBarry Smith    Calling sequence of func:
3642cdcf91faSSean Farley $    flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx);
3643325fc9f4SBarry Smith 
3644325fc9f4SBarry Smith 
3645325fc9f4SBarry Smith    Level: developer
3646325fc9f4SBarry Smith 
3647325fc9f4SBarry Smith .keywords: TS, nonlinear, set, function
3648325fc9f4SBarry Smith 
3649325fc9f4SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3650325fc9f4SBarry Smith */
3651cdcf91faSSean Farley PetscErrorCode  TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx)
3652325fc9f4SBarry Smith {
3653325fc9f4SBarry Smith   PetscErrorCode    ierr;
3654325fc9f4SBarry Smith   TSMatlabContext *sctx;
3655325fc9f4SBarry Smith 
3656325fc9f4SBarry Smith   PetscFunctionBegin;
3657325fc9f4SBarry Smith   /* currently sctx is memory bleed */
3658325fc9f4SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3659325fc9f4SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3660325fc9f4SBarry Smith   /*
3661325fc9f4SBarry Smith      This should work, but it doesn't
3662325fc9f4SBarry Smith   sctx->ctx = ctx;
3663325fc9f4SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
3664325fc9f4SBarry Smith   */
3665325fc9f4SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
3666cdcf91faSSean Farley   ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
3667325fc9f4SBarry Smith   PetscFunctionReturn(0);
3668325fc9f4SBarry Smith }
3669325fc9f4SBarry Smith 
3670b5b1a830SBarry Smith #undef __FUNCT__
3671b5b1a830SBarry Smith #define __FUNCT__ "TSMonitor_Matlab"
3672b5b1a830SBarry Smith /*
3673b5b1a830SBarry Smith    TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab().
3674b5b1a830SBarry Smith 
3675b5b1a830SBarry Smith    Collective on TS
3676b5b1a830SBarry Smith 
3677b5b1a830SBarry Smith .seealso: TSSetFunction(), TSGetFunction()
3678b5b1a830SBarry Smith @*/
3679cdcf91faSSean Farley PetscErrorCode  TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx)
3680b5b1a830SBarry Smith {
3681b5b1a830SBarry Smith   PetscErrorCode  ierr;
3682b5b1a830SBarry Smith   TSMatlabContext *sctx = (TSMatlabContext *)ctx;
3683a530c242SBarry Smith   int             nlhs = 1,nrhs = 6;
3684b5b1a830SBarry Smith   mxArray         *plhs[1],*prhs[6];
3685b5b1a830SBarry Smith   long long int   lx = 0,ls = 0;
3686b5b1a830SBarry Smith 
3687b5b1a830SBarry Smith   PetscFunctionBegin;
3688cdcf91faSSean Farley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3689b5b1a830SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,4);
3690b5b1a830SBarry Smith 
3691cdcf91faSSean Farley   ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr);
3692b5b1a830SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3693b5b1a830SBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
3694b5b1a830SBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)it);
3695b5b1a830SBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)time);
3696b5b1a830SBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lx);
3697b5b1a830SBarry Smith   prhs[4] =  mxCreateString(sctx->funcname);
3698b5b1a830SBarry Smith   prhs[5] =  sctx->ctx;
3699b5b1a830SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr);
3700b5b1a830SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3701b5b1a830SBarry Smith   mxDestroyArray(prhs[0]);
3702b5b1a830SBarry Smith   mxDestroyArray(prhs[1]);
3703b5b1a830SBarry Smith   mxDestroyArray(prhs[2]);
3704b5b1a830SBarry Smith   mxDestroyArray(prhs[3]);
3705b5b1a830SBarry Smith   mxDestroyArray(prhs[4]);
3706b5b1a830SBarry Smith   mxDestroyArray(plhs[0]);
3707b5b1a830SBarry Smith   PetscFunctionReturn(0);
3708b5b1a830SBarry Smith }
3709b5b1a830SBarry Smith 
3710b5b1a830SBarry Smith 
3711b5b1a830SBarry Smith #undef __FUNCT__
3712b5b1a830SBarry Smith #define __FUNCT__ "TSMonitorSetMatlab"
3713b5b1a830SBarry Smith /*
3714b5b1a830SBarry Smith    TSMonitorSetMatlab - Sets the monitor function from Matlab
3715b5b1a830SBarry Smith 
3716b5b1a830SBarry Smith    Level: developer
3717b5b1a830SBarry Smith 
3718b5b1a830SBarry Smith .keywords: TS, nonlinear, set, function
3719b5b1a830SBarry Smith 
3720b5b1a830SBarry Smith .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction()
3721b5b1a830SBarry Smith */
3722cdcf91faSSean Farley PetscErrorCode  TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx)
3723b5b1a830SBarry Smith {
3724b5b1a830SBarry Smith   PetscErrorCode    ierr;
3725b5b1a830SBarry Smith   TSMatlabContext *sctx;
3726b5b1a830SBarry Smith 
3727b5b1a830SBarry Smith   PetscFunctionBegin;
3728b5b1a830SBarry Smith   /* currently sctx is memory bleed */
3729b5b1a830SBarry Smith   ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr);
3730b5b1a830SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3731b5b1a830SBarry Smith   /*
3732b5b1a830SBarry Smith      This should work, but it doesn't
3733b5b1a830SBarry Smith   sctx->ctx = ctx;
3734b5b1a830SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
3735b5b1a830SBarry Smith   */
3736b5b1a830SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
3737cdcf91faSSean Farley   ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
3738b5b1a830SBarry Smith   PetscFunctionReturn(0);
3739b5b1a830SBarry Smith }
3740325fc9f4SBarry Smith #endif
3741