mirror of
https://github.com/slsdetectorgroup/aare.git
synced 2025-06-12 15:27:13 +02:00
Added fitting with lmfit (#128)
- added stand alone fitting using: https://jugit.fz-juelich.de/mlz/lmfit.git - fit_gaus, fit_pol1 with and without errors - multi threaded fitting --------- Co-authored-by: JulianHeymes <julian.heymes@psi.ch>
This commit is contained in:
223
python/src/fit.hpp
Normal file
223
python/src/fit.hpp
Normal file
@ -0,0 +1,223 @@
|
||||
#include <cstdint>
|
||||
#include <filesystem>
|
||||
#include <pybind11/pybind11.h>
|
||||
#include <pybind11/stl.h>
|
||||
#include <pybind11/stl_bind.h>
|
||||
|
||||
#include "aare/Fit.hpp"
|
||||
|
||||
namespace py = pybind11;
|
||||
|
||||
void define_fit_bindings(py::module &m) {
|
||||
|
||||
// TODO! Evaluate without converting to double
|
||||
m.def(
|
||||
"gaus",
|
||||
[](py::array_t<double, py::array::c_style | py::array::forcecast> x,
|
||||
py::array_t<double, py::array::c_style | py::array::forcecast> par) {
|
||||
auto x_view = make_view_1d(x);
|
||||
auto par_view = make_view_1d(par);
|
||||
auto y = new NDArray<double, 1>{aare::func::gaus(x_view, par_view)};
|
||||
return return_image_data(y);
|
||||
},
|
||||
R"(
|
||||
Evaluate a 1D Gaussian function for all points in x using parameters par.
|
||||
|
||||
Parameters
|
||||
----------
|
||||
x : array_like
|
||||
The points at which to evaluate the Gaussian function.
|
||||
par : array_like
|
||||
The parameters of the Gaussian function. The first element is the amplitude, the second element is the mean, and the third element is the standard deviation.
|
||||
)", py::arg("x"), py::arg("par"));
|
||||
|
||||
m.def(
|
||||
"pol1",
|
||||
[](py::array_t<double, py::array::c_style | py::array::forcecast> x,
|
||||
py::array_t<double, py::array::c_style | py::array::forcecast> par) {
|
||||
auto x_view = make_view_1d(x);
|
||||
auto par_view = make_view_1d(par);
|
||||
auto y = new NDArray<double, 1>{aare::func::pol1(x_view, par_view)};
|
||||
return return_image_data(y);
|
||||
},
|
||||
R"(
|
||||
Evaluate a 1D polynomial function for all points in x using parameters par. (p0+p1*x)
|
||||
|
||||
Parameters
|
||||
----------
|
||||
x : array_like
|
||||
The points at which to evaluate the polynomial function.
|
||||
par : array_like
|
||||
The parameters of the polynomial function. The first element is the intercept, and the second element is the slope.
|
||||
)", py::arg("x"), py::arg("par"));
|
||||
|
||||
m.def(
|
||||
"fit_gaus",
|
||||
[](py::array_t<double, py::array::c_style | py::array::forcecast> x,
|
||||
py::array_t<double, py::array::c_style | py::array::forcecast> y,
|
||||
int n_threads) {
|
||||
if (y.ndim() == 3) {
|
||||
auto par = new NDArray<double, 3>{};
|
||||
auto y_view = make_view_3d(y);
|
||||
auto x_view = make_view_1d(x);
|
||||
*par = aare::fit_gaus(x_view, y_view, n_threads);
|
||||
return return_image_data(par);
|
||||
} else if (y.ndim() == 1) {
|
||||
auto par = new NDArray<double, 1>{};
|
||||
auto y_view = make_view_1d(y);
|
||||
auto x_view = make_view_1d(x);
|
||||
*par = aare::fit_gaus(x_view, y_view);
|
||||
return return_image_data(par);
|
||||
} else {
|
||||
throw std::runtime_error("Data must be 1D or 3D");
|
||||
}
|
||||
},
|
||||
R"(
|
||||
Fit a 1D Gaussian to data.
|
||||
|
||||
Parameters
|
||||
----------
|
||||
x : array_like
|
||||
The x values.
|
||||
y : array_like
|
||||
The y values.
|
||||
n_threads : int, optional
|
||||
The number of threads to use. Default is 4.
|
||||
)",
|
||||
py::arg("x"), py::arg("y"), py::arg("n_threads") = 4);
|
||||
|
||||
m.def(
|
||||
"fit_gaus",
|
||||
[](py::array_t<double, py::array::c_style | py::array::forcecast> x,
|
||||
py::array_t<double, py::array::c_style | py::array::forcecast> y,
|
||||
py::array_t<double, py::array::c_style | py::array::forcecast>
|
||||
y_err, int n_threads) {
|
||||
if (y.ndim() == 3) {
|
||||
// Allocate memory for the output
|
||||
// Need to have pointers to allow python to manage
|
||||
// the memory
|
||||
auto par = new NDArray<double, 3>({y.shape(0), y.shape(1), 3});
|
||||
auto par_err =
|
||||
new NDArray<double, 3>({y.shape(0), y.shape(1), 3});
|
||||
|
||||
auto y_view = make_view_3d(y);
|
||||
auto y_view_err = make_view_3d(y_err);
|
||||
auto x_view = make_view_1d(x);
|
||||
aare::fit_gaus(x_view, y_view, y_view_err, par->view(),
|
||||
par_err->view(), n_threads);
|
||||
// return return_image_data(par);
|
||||
return py::make_tuple(return_image_data(par),
|
||||
return_image_data(par_err));
|
||||
} else if (y.ndim() == 1) {
|
||||
// Allocate memory for the output
|
||||
// Need to have pointers to allow python to manage
|
||||
// the memory
|
||||
auto par = new NDArray<double, 1>({3});
|
||||
auto par_err = new NDArray<double, 1>({3});
|
||||
|
||||
// Decode the numpy arrays
|
||||
auto y_view = make_view_1d(y);
|
||||
auto y_view_err = make_view_1d(y_err);
|
||||
auto x_view = make_view_1d(x);
|
||||
|
||||
aare::fit_gaus(x_view, y_view, y_view_err, par->view(),
|
||||
par_err->view());
|
||||
return py::make_tuple(return_image_data(par),
|
||||
return_image_data(par_err));
|
||||
} else {
|
||||
throw std::runtime_error("Data must be 1D or 3D");
|
||||
}
|
||||
},
|
||||
R"(
|
||||
Fit a 1D Gaussian to data with error estimates.
|
||||
|
||||
Parameters
|
||||
----------
|
||||
x : array_like
|
||||
The x values.
|
||||
y : array_like
|
||||
The y values.
|
||||
y_err : array_like
|
||||
The error in the y values.
|
||||
n_threads : int, optional
|
||||
The number of threads to use. Default is 4.
|
||||
)",
|
||||
py::arg("x"), py::arg("y"), py::arg("y_err"), py::arg("n_threads") = 4);
|
||||
|
||||
m.def(
|
||||
"fit_pol1",
|
||||
[](py::array_t<double, py::array::c_style | py::array::forcecast> x,
|
||||
py::array_t<double, py::array::c_style | py::array::forcecast> y,
|
||||
int n_threads) {
|
||||
if (y.ndim() == 3) {
|
||||
auto par = new NDArray<double, 3>{};
|
||||
|
||||
auto x_view = make_view_1d(x);
|
||||
auto y_view = make_view_3d(y);
|
||||
*par = aare::fit_pol1(x_view, y_view, n_threads);
|
||||
return return_image_data(par);
|
||||
} else if (y.ndim() == 1) {
|
||||
auto par = new NDArray<double, 1>{};
|
||||
auto x_view = make_view_1d(x);
|
||||
auto y_view = make_view_1d(y);
|
||||
*par = aare::fit_pol1(x_view, y_view);
|
||||
return return_image_data(par);
|
||||
} else {
|
||||
throw std::runtime_error("Data must be 1D or 3D");
|
||||
}
|
||||
},
|
||||
py::arg("x"), py::arg("y"), py::arg("n_threads") = 4);
|
||||
|
||||
m.def(
|
||||
"fit_pol1",
|
||||
[](py::array_t<double, py::array::c_style | py::array::forcecast> x,
|
||||
py::array_t<double, py::array::c_style | py::array::forcecast> y,
|
||||
py::array_t<double, py::array::c_style | py::array::forcecast>
|
||||
y_err, int n_threads) {
|
||||
if (y.ndim() == 3) {
|
||||
auto par =
|
||||
new NDArray<double, 3>({y.shape(0), y.shape(1), 2});
|
||||
auto par_err =
|
||||
new NDArray<double, 3>({y.shape(0), y.shape(1), 2});
|
||||
|
||||
auto y_view = make_view_3d(y);
|
||||
auto y_view_err = make_view_3d(y_err);
|
||||
auto x_view = make_view_1d(x);
|
||||
|
||||
aare::fit_pol1(x_view, y_view,y_view_err, par->view(),
|
||||
par_err->view(), n_threads);
|
||||
return py::make_tuple(return_image_data(par),
|
||||
return_image_data(par_err));
|
||||
|
||||
} else if (y.ndim() == 1) {
|
||||
auto par = new NDArray<double, 1>({2});
|
||||
auto par_err = new NDArray<double, 1>({2});
|
||||
|
||||
auto y_view = make_view_1d(y);
|
||||
auto y_view_err = make_view_1d(y_err);
|
||||
auto x_view = make_view_1d(x);
|
||||
|
||||
aare::fit_pol1(x_view, y_view, y_view_err, par->view(),
|
||||
par_err->view());
|
||||
return py::make_tuple(return_image_data(par),
|
||||
return_image_data(par_err));
|
||||
} else {
|
||||
throw std::runtime_error("Data must be 1D or 3D");
|
||||
}
|
||||
},
|
||||
R"(
|
||||
Fit a 1D polynomial to data with error estimates.
|
||||
|
||||
Parameters
|
||||
----------
|
||||
x : array_like
|
||||
The x values.
|
||||
y : array_like
|
||||
The y values.
|
||||
y_err : array_like
|
||||
The error in the y values.
|
||||
n_threads : int, optional
|
||||
The number of threads to use. Default is 4.
|
||||
)",
|
||||
py::arg("x"), py::arg("y"), py::arg("y_err"), py::arg("n_threads") = 4);
|
||||
}
|
@ -8,6 +8,7 @@
|
||||
#include "pedestal.hpp"
|
||||
#include "cluster.hpp"
|
||||
#include "cluster_file.hpp"
|
||||
#include "fit.hpp"
|
||||
|
||||
//Pybind stuff
|
||||
#include <pybind11/pybind11.h>
|
||||
@ -29,5 +30,6 @@ PYBIND11_MODULE(_aare, m) {
|
||||
define_cluster_file_io_bindings(m);
|
||||
define_cluster_collector_bindings(m);
|
||||
define_cluster_file_sink_bindings(m);
|
||||
define_fit_bindings(m);
|
||||
|
||||
}
|
@ -39,65 +39,6 @@ template <typename T> py::array return_vector(std::vector<T> *vec) {
|
||||
free_when_done); // numpy array references this parent
|
||||
}
|
||||
|
||||
// template <typename Reader> py::array do_read(Reader &r, size_t n_frames) {
|
||||
// py::array image;
|
||||
// if (n_frames == 0)
|
||||
// n_frames = r.total_frames();
|
||||
|
||||
// std::array<ssize_t, 3> shape{static_cast<ssize_t>(n_frames), r.rows(),
|
||||
// r.cols()};
|
||||
// const uint8_t item_size = r.bytes_per_pixel();
|
||||
// if (item_size == 1) {
|
||||
// image = py::array_t<uint8_t, py::array::c_style | py::array::forcecast>(
|
||||
// shape);
|
||||
// } else if (item_size == 2) {
|
||||
// image =
|
||||
// py::array_t<uint16_t, py::array::c_style | py::array::forcecast>(
|
||||
// shape);
|
||||
// } else if (item_size == 4) {
|
||||
// image =
|
||||
// py::array_t<uint32_t, py::array::c_style | py::array::forcecast>(
|
||||
// shape);
|
||||
// }
|
||||
// r.read_into(reinterpret_cast<std::byte *>(image.mutable_data()), n_frames);
|
||||
// return image;
|
||||
// }
|
||||
|
||||
// py::array return_frame(pl::Frame *ptr) {
|
||||
// py::capsule free_when_done(ptr, [](void *f) {
|
||||
// pl::Frame *foo = reinterpret_cast<pl::Frame *>(f);
|
||||
// delete foo;
|
||||
// });
|
||||
|
||||
// const uint8_t item_size = ptr->bytes_per_pixel();
|
||||
// std::vector<ssize_t> shape;
|
||||
// for (auto val : ptr->shape())
|
||||
// if (val > 1)
|
||||
// shape.push_back(val);
|
||||
|
||||
// std::vector<ssize_t> strides;
|
||||
// if (shape.size() == 1)
|
||||
// strides.push_back(item_size);
|
||||
// else if (shape.size() == 2) {
|
||||
// strides.push_back(item_size * shape[1]);
|
||||
// strides.push_back(item_size);
|
||||
// }
|
||||
|
||||
// if (item_size == 1)
|
||||
// return py::array_t<uint8_t>(
|
||||
// shape, strides,
|
||||
// reinterpret_cast<uint8_t *>(ptr->data()), free_when_done);
|
||||
// else if (item_size == 2)
|
||||
// return py::array_t<uint16_t>(shape, strides,
|
||||
// reinterpret_cast<uint16_t *>(ptr->data()),
|
||||
// free_when_done);
|
||||
// else if (item_size == 4)
|
||||
// return py::array_t<uint32_t>(shape, strides,
|
||||
// reinterpret_cast<uint32_t *>(ptr->data()),
|
||||
// free_when_done);
|
||||
// return {};
|
||||
// }
|
||||
|
||||
// todo rewrite generic
|
||||
template <class T, int Flags> auto get_shape_3d(py::array_t<T, Flags> arr) {
|
||||
return aare::Shape<3>{arr.shape(0), arr.shape(1), arr.shape(2)};
|
||||
@ -111,6 +52,13 @@ template <class T, int Flags> auto get_shape_2d(py::array_t<T, Flags> arr) {
|
||||
return aare::Shape<2>{arr.shape(0), arr.shape(1)};
|
||||
}
|
||||
|
||||
template <class T, int Flags> auto get_shape_1d(py::array_t<T, Flags> arr) {
|
||||
return aare::Shape<1>{arr.shape(0)};
|
||||
}
|
||||
|
||||
template <class T, int Flags> auto make_view_2d(py::array_t<T, Flags> arr) {
|
||||
return aare::NDView<T, 2>(arr.mutable_data(), get_shape_2d<T, Flags>(arr));
|
||||
}
|
||||
template <class T, int Flags> auto make_view_1d(py::array_t<T, Flags> arr) {
|
||||
return aare::NDView<T, 1>(arr.mutable_data(), get_shape_1d<T, Flags>(arr));
|
||||
}
|
Reference in New Issue
Block a user