Skip to content

Commit

Permalink
pair_t -> scale_t
Browse files Browse the repository at this point in the history
Rename a few structs to clarify the semantic meanings. E.g.
`scale_t<Micron>` implies voxel sizes in micrometer. Similarly,
`scale_t<Micron>.xy` implies a square pixel.
  • Loading branch information
antonysigma committed Nov 18, 2024
1 parent 1fd6ce8 commit a16679e
Show file tree
Hide file tree
Showing 2 changed files with 13 additions and 13 deletions.
6 changes: 3 additions & 3 deletions microsc-psf/inc/make_psf.h
Original file line number Diff line number Diff line change
Expand Up @@ -43,11 +43,11 @@ struct precision_li2017_t {
};

template <typename T>
struct pair_t {
T x{};
struct scale_t {
T xy{};
T z{};
};

arma::Cube<double> makePSF(microscope_params_t, pair_t<Micron> voxel, pair_t<uint32_t> volume,
arma::Cube<double> makePSF(microscope_params_t, scale_t<Micron> voxel, scale_t<uint32_t> volume,
Micron wavelength = 0.530_um, precision_li2017_t = {});
} // namespace microsc_psf
20 changes: 10 additions & 10 deletions microsc-psf/src/main.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -23,27 +23,27 @@ iota(uint32_t N) {
/** Piecewise linear interpolation. */
arma::Cube<double>
cylToRectTransform(const arma::mat& PSF0, const arma::vec& R,
microsc_psf::pair_t<uint32_t> volume) {
microsc_psf::scale_t<uint32_t> volume) {
using namespace arma;
assert(PSF0.n_cols == uint32_t(volume.z));

vec rPixel;
{
const double x0 = (volume.x - 1) / 2.0;
const double x0 = (volume.xy - 1) / 2.0;

// meshgrid(0:nx - x0, 0:ny - y0)
mat X(volume.x, volume.x);
X.each_row() = (iota(volume.x) - x0).t();
mat X(volume.xy, volume.xy);
X.each_row() = (iota(volume.xy) - x0).t();
#define Y X.t()

rPixel = sqrt(X % X + Y % Y).as_col();
}

Cube<double> PSF(volume.x, volume.x, volume.z);
Cube<double> PSF(volume.xy, volume.xy, volume.z);
#pragma omp parallel for
for (uint32_t zi = 0; zi < volume.z; zi++) {
// Memory map the PSF slice to an 1D vector without memory copy.
vec interpolated{&(PSF(0, 0, zi)), static_cast<uint64_t>(volume.x) * volume.x, false, true};
vec interpolated{&(PSF(0, 0, zi)), static_cast<uint64_t>(volume.xy) * volume.xy, false, true};

interp1(R, PSF0.col(zi), rPixel, interpolated, "linear", 0.0);
}
Expand All @@ -55,18 +55,18 @@ cylToRectTransform(const arma::mat& PSF0, const arma::vec& R,
namespace microsc_psf {

arma::Cube<double>
makePSF(microscope_params_t params, pair_t<Micron> voxel, pair_t<uint32_t> volume,
makePSF(microscope_params_t params, scale_t<Micron> voxel, scale_t<uint32_t> volume,
Micron wavelength, precision_li2017_t precision) {
using ::units::literals::operator""_m;

const double x0 = (volume.x - 1) / 2.0;
const double x0 = (volume.xy - 1) / 2.0;
const double y0 = x0;
const double z0 = (volume.z - 1) / 2.0;

using namespace arma;

// Max radius is the length of the diagonal of the volume xy plane.
const double max_radius = round(abs(cx_double{volume.x - x0, volume.x - y0})) + 1;
const double max_radius = round(abs(cx_double{volume.xy - x0, volume.xy - y0})) + 1;

const float max_rho = std::min({
params.NA, //
Expand Down Expand Up @@ -115,7 +115,7 @@ makePSF(microscope_params_t params, pair_t<Micron> voxel, pair_t<uint32_t> volum
// See equation 5 in Li, Xue, and Blu 2017.
mat Ele;
{
const vec A = k0 * params.NA * Meter(voxel.x) * R;
const vec A = k0 * params.NA * Meter(voxel.xy) * R;

// bsxfun(@minus, an2', A2);
mat domin(A.n_elem, scaling_factor.n_elem);
Expand Down

0 comments on commit a16679e

Please sign in to comment.