Skip to content
Snippets Groups Projects
Commit 9e4f87fb authored by akuesters's avatar akuesters Committed by Benjamin Cumming
Browse files

Python wrapper: complete schedule by adding events(t0,t1) (#826)

* add interface for Python users to query schedule event times
* add schedule event tests and assertions
* add schedule events to doc
parent 51f83b42
No related branches found
No related tags found
No related merge requests found
......@@ -193,6 +193,10 @@ Event Generator and Schedules
No events delivered after this time [ms].
Must be non-negative or ``None``.
.. function:: events(t0, t1)
Returns a view of monotonically increasing time values in the half-open interval [t0, t1).
.. class:: explicit_schedule
Describes an explicit schedule at a predetermined (sorted) sequence of :attr:`times`.
......@@ -207,6 +211,10 @@ Event Generator and Schedules
The list of non-negative times [ms].
.. function:: events(t0, t1)
Returns a view of monotonically increasing time values in the half-open interval [t0, t1).
.. class:: poisson_schedule
Describes a schedule according to a Poisson process.
......@@ -230,6 +238,10 @@ Event Generator and Schedules
The seed for the random number generator.
.. function:: events(t0, t1)
Returns a view of monotonically increasing time values in the half-open interval [t0, t1).
An example of an event generator reads as follows:
.. container:: example-code
......
......@@ -3,6 +3,7 @@
#include <arbor/util/optional.hpp>
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
#include "conversion.hpp"
#include "schedule.hpp"
......@@ -28,6 +29,10 @@ std::ostream& operator<<(std::ostream& o, const poisson_schedule_shim& p) {
<< ", seed " << p.seed << ">";
};
static std::vector<arb::time_type> as_vector(std::pair<const arb::time_type*, const arb::time_type*> ts) {
return std::vector<arb::time_type>(ts.first, ts.second);
}
//
// regular_schedule shim
//
......@@ -76,6 +81,15 @@ arb::schedule regular_schedule_shim::schedule() const {
tstop.value_or(arb::terminal_time));
}
std::vector<arb::time_type> regular_schedule_shim::events(arb::time_type t0, arb::time_type t1) {
pyarb::assert_throw(is_nonneg()(t0), "t0 must be a non-negative number");
pyarb::assert_throw(is_nonneg()(t1), "t1 must be a non-negative number");
arb::schedule sched = regular_schedule_shim::schedule();
return as_vector(sched.events(t0, t1));
}
//
// explicit_schedule shim
//
......@@ -109,6 +123,15 @@ arb::schedule explicit_schedule_shim::schedule() const {
return arb::explicit_schedule(times);
}
std::vector<arb::time_type> explicit_schedule_shim::events(arb::time_type t0, arb::time_type t1) {
pyarb::assert_throw(is_nonneg()(t0), "t0 must be a non-negative number");
pyarb::assert_throw(is_nonneg()(t1), "t1 must be a non-negative number");
arb::schedule sched = explicit_schedule_shim::schedule();
return as_vector(sched.events(t0, t1));
}
//
// poisson_schedule shim
//
......@@ -146,6 +169,15 @@ arb::schedule poisson_schedule_shim::schedule() const {
return arb::poisson_schedule(tstart, freq/1000., rng_type(seed));
}
std::vector<arb::time_type> poisson_schedule_shim::events(arb::time_type t0, arb::time_type t1) {
pyarb::assert_throw(is_nonneg()(t0), "t0 must be a non-negative number");
pyarb::assert_throw(is_nonneg()(t1), "t1 must be a non-negative number");
arb::schedule sched = poisson_schedule_shim::schedule();
return as_vector(sched.events(t0, t1));
}
void register_schedules(pybind11::module& m) {
using namespace pybind11::literals;
using time_type = arb::time_type;
......@@ -167,6 +199,8 @@ void register_schedules(pybind11::module& m) {
"No events delivered after this time [ms].")
.def_property("dt", &regular_schedule_shim::get_dt, &regular_schedule_shim::set_dt,
"The interval between time points [ms].")
.def("events", &regular_schedule_shim::events,
"A view of monotonically increasing time values in the half-open interval [t0, t1).")
.def("__str__", util::to_string<regular_schedule_shim>)
.def("__repr__", util::to_string<regular_schedule_shim>);
......@@ -183,6 +217,8 @@ void register_schedules(pybind11::module& m) {
" times: A list of times [ms], [] by default.")
.def_property("times", &explicit_schedule_shim::get_times, &explicit_schedule_shim::set_times,
"A list of times [ms].")
.def("events", &explicit_schedule_shim::events,
"A view of monotonically increasing time values in the half-open interval [t0, t1).")
.def("__str__", util::to_string<explicit_schedule_shim>)
.def("__repr__", util::to_string<explicit_schedule_shim>);
......@@ -203,6 +239,8 @@ void register_schedules(pybind11::module& m) {
"The expected frequency [Hz].")
.def_readwrite("seed", &poisson_schedule_shim::seed,
"The seed for the random number generator.")
.def("events", &poisson_schedule_shim::events,
"A view of monotonically increasing time values in the half-open interval [t0, t1).")
.def("__str__", util::to_string<poisson_schedule_shim>)
.def("__repr__", util::to_string<poisson_schedule_shim>);
}
......
#pragma once
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
#include <arbor/schedule.hpp>
#include <arbor/common_types.hpp>
......@@ -34,6 +35,8 @@ struct regular_schedule_shim {
opt_time_type get_tstop() const;
arb::schedule schedule() const;
std::vector<arb::time_type> events(arb::time_type t0, arb::time_type t1);
};
// A Python shim for arb::explicit_schedule.
......@@ -51,6 +54,8 @@ struct explicit_schedule_shim {
std::vector<arb::time_type> get_times() const;
arb::schedule schedule() const;
std::vector<arb::time_type> events(arb::time_type t0, arb::time_type t1);
};
// A Python shim for arb::poisson_schedule.
......@@ -74,6 +79,8 @@ struct poisson_schedule_shim {
arb::time_type get_freq() const;
arb::schedule schedule() const;
std::vector<arb::time_type> events(arb::time_type t0, arb::time_type t1);
};
}
......@@ -35,9 +35,16 @@ class RegularSchedule(unittest.TestCase):
rs.dt = 0.5
rs.tstop = 42.
self.assertEqual(rs.tstart, 17.)
self.assertAlmostEqual(rs.dt, 0.5)
self.assertAlmostEqual(rs.dt, 0.5, places = 1)
self.assertEqual(rs.tstop, 42.)
def test_events_regular_schedule(self):
expected = [0, 0.25, 0.5, 0.75, 1.0]
rs = arb.regular_schedule(tstart = 0., dt = 0.25, tstop = 1.25)
self.assertEqual(expected, rs.events(0., 1.25))
self.assertEqual(expected, rs.events(0., 5.))
self.assertEqual([], rs.events(5., 10.))
def test_exceptions_regular_schedule(self):
with self.assertRaisesRegex(RuntimeError,
"tstart must be a non-negative number, or None"):
......@@ -52,6 +59,14 @@ class RegularSchedule(unittest.TestCase):
with self.assertRaisesRegex(RuntimeError,
"tstop must be a non-negative number, or None"):
arb.regular_schedule(tstop = 'tstop')
with self.assertRaisesRegex(RuntimeError,
"t0 must be a non-negative number"):
rs = arb.regular_schedule(0., 1., 10.)
rs.events(-1, 0)
with self.assertRaisesRegex(RuntimeError,
"t1 must be a non-negative number"):
rs = arb.regular_schedule(0., 1., 10.)
rs.events(0, -10)
class ExplicitSchedule(unittest.TestCase):
def test_times_contor_explicit_schedule(self):
......@@ -63,6 +78,16 @@ class ExplicitSchedule(unittest.TestCase):
es.times = [42, 43, 44, 55.5, 100]
self.assertEqual(es.times, [42, 43, 44, 55.5, 100])
def test_events_explicit_schedule(self):
times = [0.1, 0.3, 1.0, 2.2, 1.25, 1.7]
expected = [0.1, 0.3, 1.0]
es = arb.explicit_schedule(times)
for i in range(len(expected)):
self.assertAlmostEqual(expected[i], es.events(0., 1.25)[i], places = 2)
expected = [0.3, 1.0, 1.25, 1.7]
for i in range(len(expected)):
self.assertAlmostEqual(expected[i], es.events(0.3, 1.71)[i], places = 2)
def test_exceptions_explicit_schedule(self):
with self.assertRaisesRegex(RuntimeError,
"explicit time schedule cannot contain negative values"):
......@@ -73,6 +98,14 @@ class ExplicitSchedule(unittest.TestCase):
arb.explicit_schedule([None])
with self.assertRaises(TypeError):
arb.explicit_schedule([[1,2,3]])
with self.assertRaisesRegex(RuntimeError,
"t0 must be a non-negative number"):
rs = arb.regular_schedule()
rs.events(-1., 1.)
with self.assertRaisesRegex(RuntimeError,
"t1 must be a non-negative number"):
rs = arb.regular_schedule()
rs.events(1., -1.)
class PoissonSchedule(unittest.TestCase):
def test_freq_seed_contor_poisson_schedule(self):
......@@ -95,6 +128,15 @@ class PoissonSchedule(unittest.TestCase):
self.assertAlmostEqual(ps.freq, 5.5)
self.assertEqual(ps.seed, 83)
def test_events_poisson_schedule(self):
expected = [17.4107, 502.074, 506.111, 597.116]
ps = arb.poisson_schedule(0., 10., 0)
for i in range(len(expected)):
self.assertAlmostEqual(expected[i], ps.events(0., 600.)[i], places = 3)
expected = [5030.22, 5045.75, 5069.84, 5091.56, 5182.17, 5367.3, 5566.73, 5642.13, 5719.85, 5796, 5808.33]
for i in range(len(expected)):
self.assertAlmostEqual(expected[i], ps.events(5000., 6000.)[i], places = 2)
def test_exceptions_poisson_schedule(self):
with self.assertRaisesRegex(RuntimeError,
"tstart must be a non-negative number"):
......@@ -116,6 +158,14 @@ class PoissonSchedule(unittest.TestCase):
arb.poisson_schedule(seed = 'seed')
with self.assertRaises(TypeError):
arb.poisson_schedule(seed = None)
with self.assertRaisesRegex(RuntimeError,
"t0 must be a non-negative number"):
ps = arb.poisson_schedule()
ps.events(-1., 1.)
with self.assertRaisesRegex(RuntimeError,
"t1 must be a non-negative number"):
ps = arb.poisson_schedule()
ps.events(1., -1.)
def suite():
# specify class and test functions in tuple (here: all tests starting with 'test' from classes RegularSchedule, ExplicitSchedule and PoissonSchedule
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment