Skip to content

Commit

Permalink
update to pybind11 v2.2.0
Browse files Browse the repository at this point in the history
  • Loading branch information
michaelbuehlmann committed Aug 31, 2017
1 parent cdcf1a3 commit 4b62406
Show file tree
Hide file tree
Showing 6 changed files with 16 additions and 22 deletions.
2 changes: 1 addition & 1 deletion external/pybind11
Submodule pybind11 updated 160 files
10 changes: 4 additions & 6 deletions python/py_basictypes.cpp
Expand Up @@ -16,9 +16,9 @@ namespace py = pybind11;
using namespace catana;


PYBIND11_PLUGIN(basictypes)
PYBIND11_MODULE(basictypes, m)
{
py::module m("basictypes", "python bindings of basic types used in CatAna");
m.doc() = "python bindings of basic types used in CatAna";

// Binding for the "Point" class, with cartesian and spherical output methods
py::class_<Point>(m, "Point", py::buffer_protocol(), R"pbdoc(
Expand Down Expand Up @@ -125,11 +125,11 @@ Create a PointContainer from a (numpy) array containing the coordinates
// Buffer definition, so that we can call numpy.array(point_container) in python
.def_buffer([](PointContainer& oc) -> py::buffer_info {
return py::buffer_info(
(double*) oc.data(),
&(oc.data()->r),
sizeof(double),
py::format_descriptor<double>::format(),
2,
{ oc.size(), 3},
std::array<long int, 2>{static_cast<long int>(oc.size()), 3},
{ sizeof(double) * 3, sizeof(double)}
);
})
Expand Down Expand Up @@ -215,6 +215,4 @@ Evaluate the FunctionInterpolator at the given values x.
Parameters:
x (numpy.ndarray[float]): evaluation points
)pbdoc", py::arg("x"));

return m.ptr();
}
6 changes: 2 additions & 4 deletions python/py_besseltools.cpp
Expand Up @@ -9,8 +9,8 @@ namespace py = pybind11;
using namespace catana;


PYBIND11_PLUGIN(besseltools) {
py::module m("besseltools", "python binding to besseltools library (part of CatAna)");
PYBIND11_MODULE(besseltools, m) {
m.doc() = "python binding to besseltools library (part of CatAna)";

py::class_<besseltools::SphericalBesselZeros>(m, "SphericalBesselZeros", R"pbdoc(
A root generator for spherical Bessel functions.
Expand Down Expand Up @@ -82,6 +82,4 @@ Computes :math:`\int_0^{r_{max}} \; w(r) \; j_l(k_1 r) \; j_l(k_2 r) \;\; dr`
k2 (float): frequency of other spherical Bessel function
)pbdoc", py::arg("w"), py::arg("l"), py::arg("r_max"), py::arg("k1"), py::arg("k2")
);

return m.ptr();
}
6 changes: 2 additions & 4 deletions python/py_decomposition.cpp
Expand Up @@ -10,8 +10,8 @@
namespace py = pybind11;
using namespace catana;

PYBIND11_PLUGIN(decomposition) {
py::module m("decomposition", "python binding for SFB decomposition of particle positions (part of CatAna)");
PYBIND11_MODULE(decomposition, m) {
m.doc() = "python binding for SFB decomposition of particle positions (part of CatAna)";

// Initialization of random numbers. Once with given seed, once random seed
m.def("init_random", py::overload_cast<unsigned int>(&init_random), R"pbdoc(
Expand Down Expand Up @@ -151,6 +151,4 @@ This method makes use of the fast pixelized decomposition scheme.
f_lm(k) components
)pbdoc", py::arg("pixelized_point_container"), py::arg("lmax"), py::arg("nmax"), py::arg("rmax"),
py::arg("store_flmn") = false, py::arg("verbose") = true);

return m.ptr();
}
12 changes: 5 additions & 7 deletions python/py_io.cpp
Expand Up @@ -16,8 +16,8 @@ using record_sd = io::SphericalRecord<double>;
using record_sd_3dex = io::SphericalRecord<double, io::SphericalTextFormat::THREEDEX>;


PYBIND11_PLUGIN(io) {
py::module m("io", "python binding for in/output of particle positions (part of CatAna)");
PYBIND11_MODULE(io, m) {
m.doc() = "python binding for in/output of particle positions (part of CatAna)";

// Initialization of random numbers. Once with given seed, once random seed
m.def("init_random", py::overload_cast<unsigned int> (&init_random), R"pbdoc(
Expand Down Expand Up @@ -343,8 +343,8 @@ Constructor from a HEALPix map array
sink (Sink): to where the filtered points are stored
buffer_size (int): size of the buffer (in number of points). Faster if larger, but consumes more memroy.
verbose (bool): if ``True``, print additional logging information to stdout
)pbdoc", py::arg("source"), py::arg("sink"), py::arg("buffer_size")=1000000, py::arg("verbose")=true,
py::keep_alive<1, 2>(), py::keep_alive<1, 3>())
)pbdoc", py::arg("source"), py::arg("sink"), py::arg("buffer_size")=1000000, py::arg("verbose")=true)
// py::keep_alive<1, 2>(), py::keep_alive<1, 3>()) TODO: add keep_alive back in (removed because of SIGSEGV)
.def("add_filter", &io::FilterStream::add_filter, R"pbdoc(
Add a Filter to the pipeline.
Expand All @@ -361,7 +361,7 @@ Constructor from a HEALPix map array
Parameters:
source (Source): a source instance
)pbdoc", py::arg("source"), py::keep_alive<1, 2>())
)pbdoc") //, py::arg("source"), py::keep_alive<1, 2>()) // TODO: add keep_alive back in
.def("run", &io::FilterStream::run, R"pbdoc(
Run the pipeline
Expand Down Expand Up @@ -402,6 +402,4 @@ Constructor from a HEALPix map array
Returns:
int: number of points written to the sink
)pbdoc", py::arg("temp_filename")="tmp.bin", py::arg("subsample_size")=0, py::arg("remove_temp")=true);

return m.ptr();
};
2 changes: 2 additions & 0 deletions test/python/io_test.py
Expand Up @@ -3,10 +3,12 @@
import pytest
import os


@pytest.fixture
def datadir():
return os.path.join(os.path.dirname(__file__), '..', 'test_data')


class TestText(object):
def test_source(self, datadir):
s = catana.io.CartesianTextSource(os.path.join(datadir, "mock_data.txt"))
Expand Down

0 comments on commit 4b62406

Please sign in to comment.