Skip to content

Commit

Permalink
[Reactor] Make argument to ReactorNet.step optional and deprecated
Browse files Browse the repository at this point in the history
The value of this argument has almost no effect on the integrator, and
frequently confuses users since the ReactorNet can end up at a time either
greater or less than the specified time. By removing this argument, the
distinction betwen step() and advance(t) becomes much more clear.
  • Loading branch information
speth committed Jun 11, 2015
1 parent 392b0f5 commit 8d4e9bf
Show file tree
Hide file tree
Showing 9 changed files with 39 additions and 31 deletions.
7 changes: 4 additions & 3 deletions include/cantera/zeroD/ReactorNet.h
Original file line number Diff line number Diff line change
Expand Up @@ -107,9 +107,10 @@ class ReactorNet : public Cantera::FuncEval
*/
void advance(doublereal time);

//! Advance the state of all reactors in time. Take a single timestep
//! toward *time*.
double step(doublereal time);
//! Advance the state of all reactors in time.
//! @deprecated The *time* argument to this function is deprecated and will
//! be removed after Cantera 2.3.
double step(doublereal time=-999);

//@}

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -31,7 +31,7 @@
# Rmax is the maximum relative reaction rate at any timestep
Rmax = np.zeros(gas.n_reactions)
while t < 0.02:
t = sim.step(1.0)
t = sim.step()
tt.append(1000 * t)
TT.append(r.T)
rnet = abs(gas.net_rates_of_progress)
Expand Down Expand Up @@ -73,7 +73,7 @@
tt = []
TT = []
while t < 0.02:
t = sim.step(1.0)
t = sim.step()
tt.append(1000 * t)
TT.append(r.T)

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,7 @@
net = ct.ReactorNet([r])
T = r.T
while T < 1900:
net.step(1.0)
net.step()
T = r.T

element = 'N'
Expand Down
2 changes: 1 addition & 1 deletion interfaces/cython/cantera/examples/reactors/combustor.py
Original file line number Diff line number Diff line change
Expand Up @@ -81,7 +81,7 @@
csvwriter = csv.writer(outfile)

while tnow < tfinal:
tnow = sim.step(tfinal)
tnow = sim.step()
tres = combustor.mass/v.mdot(tnow)
Tnow = combustor.T
if abs(Tnow - Tprev) > 1.0 or tnow-tprev > 2e-2:
Expand Down
10 changes: 7 additions & 3 deletions interfaces/cython/cantera/reactor.pyx
Original file line number Diff line number Diff line change
Expand Up @@ -764,10 +764,14 @@ cdef class ReactorNet:
"""
self.net.advance(t)

def step(self, double t):
def step(self, double t=-999):
"""
Take a single internal time step toward time *t* [s]. The time after
taking the step is returned.
Take a single internal time step. The time after taking the step is
returned.
.. deprecated:: 2.2
The argument *t* is deprecated and will be removed after
Cantera 2.3.
"""
return self.net.step(t)

Expand Down
2 changes: 1 addition & 1 deletion interfaces/cython/cantera/test/test_kinetics.py
Original file line number Diff line number Diff line change
Expand Up @@ -363,7 +363,7 @@ def test_dot_output(self):
net = ct.ReactorNet([r])
T = r.T
while T < 1900:
net.step(1.0)
net.step()
T = r.T

for element in ['N','C','H','O']:
Expand Down
26 changes: 13 additions & 13 deletions interfaces/cython/cantera/test/test_reactor.py
Original file line number Diff line number Diff line change
Expand Up @@ -77,7 +77,7 @@ def test_names(self):

def test_component_index(self):
self.make_reactors(n_reactors=1)
self.net.step(1.0)
self.net.step()

N0 = self.net.n_vars - self.gas1.n_species
for i, name in enumerate(self.gas1.species_names):
Expand Down Expand Up @@ -122,7 +122,7 @@ def test_timestepping(self):

while t < tEnd:
tPrev = t
t = self.net.step(tEnd)
t = self.net.step()
self.assertTrue(t - tPrev <= 1.0001 * dt_max)
self.assertNear(t, self.net.time)

Expand Down Expand Up @@ -160,7 +160,7 @@ def integrate(atol, rtol):
t = 0

while t < tEnd:
t = self.net.step(tEnd)
t = self.net.step()
nSteps += 1

return nSteps
Expand Down Expand Up @@ -426,7 +426,7 @@ def speciesMass(k):

t = 0
while t < 1.0:
t = self.net.step(1.0)
t = self.net.step()
p1 = self.r1.thermo.P
p2 = self.r2.thermo.P
self.assertNear(mdot(p1-p2), valve.mdot(t))
Expand All @@ -452,7 +452,7 @@ def test_pressure_controller(self):

t = 0
while t < 1.0:
t = self.net.step(1.0)
t = self.net.step()
self.assertNear(mdot(t), mfc.mdot(t))
dP = self.r1.thermo.P - outlet_reservoir.thermo.P
self.assertNear(mdot(t) + 1e-5 * dP, pc.mdot(t))
Expand Down Expand Up @@ -585,7 +585,7 @@ def integrate(self, tf):
i = 0
while t < tf:
i += 1
t = self.net.step(tf)
t = self.net.step()
times.append(t)
T.append(self.combustor.T)
return times, T
Expand Down Expand Up @@ -706,7 +706,7 @@ def test_component_index(self):
self.create_reactors(add_surf=True)
for (gas,net,iface,r) in ((self.gas1, self.net1, self.interface1, self.r1),
(self.gas2, self.net2, self.interface2, self.r2)):
net.step(1.0)
net.step()

N0 = net.n_vars - gas.n_species - iface.n_species
N1 = net.n_vars - iface.n_species
Expand Down Expand Up @@ -765,7 +765,7 @@ def test_nonreacting(self):
v0 = r.speed
self.assertNear(v0, 10 / r.density)
while t < 10.0:
t = net.step(10.0)
t = net.step()

self.assertNear(v0, r.speed)
self.assertNear(r.distance, v0 * t)
Expand All @@ -792,7 +792,7 @@ def test_reacting(self):
t1 = net.time
x1 = r.distance

t = net.step(1.0)
t = net.step()

v = (r.distance - x1) / (net.time - t1)
self.assertNear(r.speed, v, 1e-3)
Expand Down Expand Up @@ -952,7 +952,7 @@ def test_sensitivities2(self):

for T in (901, 905, 910, 950, 1500):
while r2.T < T:
net.step(1.0)
net.step()

S = net.sensitivities()

Expand Down Expand Up @@ -984,7 +984,7 @@ def setup():

def integrate(r, net):
while r.T < 910:
net.step(1.0)
net.step()
return net.sensitivities()

r1,net1 = setup()
Expand Down Expand Up @@ -1215,7 +1215,7 @@ def test_integrateWithStep(self):
tfinal = 0.25
self.data = []
while tnow < tfinal:
tnow = self.sim.step(tfinal)
tnow = self.sim.step()
self.data.append([tnow, self.combustor.T] +
list(self.combustor.thermo.X))

Expand Down Expand Up @@ -1280,7 +1280,7 @@ def test_integrateWithStep(self):
tfinal = 0.01
self.data = []
while tnow < tfinal:
tnow = self.sim.step(tfinal)
tnow = self.sim.step()
self.data.append([tnow,
self.r1.T, self.r2.T,
self.r1.thermo.P, self.r2.thermo.P,
Expand Down
11 changes: 5 additions & 6 deletions interfaces/matlab/toolbox/@ReactorNet/step.m
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
function t = step(r, tout)
% STEP Take one internal time step.
% t = step(r, tout)
% t = step(r)
% The integrator used to integrate the ODEs (CVODE) takes
% variable-size steps, chosen so that a specified error
% tolerance is maintained. At times when the solution is rapidly
Expand All @@ -19,7 +19,7 @@
% t = 0.0
% tout = 0.1
% while t < tout
% t = step(r, tout)
% t = step(r)
% ,,,
% <commands to save desired variables>
% ...
Expand All @@ -29,12 +29,11 @@
%
% :param r:
% Instance of class :mat:func:`ReactorNet`
% :param tout:
% End time that the integrator should take one internal time
% step towards. This end time may not be reached, or may be
% exceeded by the internal time step.
% :return:
% Network time after the internal time step. Units: s
%
if nargin == 1
tout = -999;
end

t = reactornetmethods(21, reactornet_hndl(r), tout);
6 changes: 5 additions & 1 deletion src/zeroD/ReactorNet.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -124,12 +124,16 @@ void ReactorNet::advance(doublereal time)

double ReactorNet::step(doublereal time)
{
if (time != -999) {
warn_deprecated("ReactorNet::step(t)", "The argument to this function"
" is deprecated and will be removed after Cantera 2.3.");
}
if (!m_init) {
initialize();
} else if (!m_integrator_init) {
reinitialize();
}
m_time = m_integ->step(time);
m_time = m_integ->step(m_time + 1.0);
updateState(m_integ->solution());
return m_time;
}
Expand Down

0 comments on commit 8d4e9bf

Please sign in to comment.