efe882f4b6
Build Packages / Unit tests (push) Failing after 1s
Build Packages / build:rpm (ubuntu2404_nocuda) (push) Successful in 25m52s
Build Packages / build:rpm (ubuntu2204_nocuda) (push) Successful in 29m5s
Build Packages / build:rpm (rocky8_nocuda) (push) Successful in 29m54s
Build Packages / build:rpm (rocky8) (push) Successful in 31m55s
Build Packages / build:rpm (rocky8_sls9) (push) Successful in 32m12s
Build Packages / build:rpm (rocky9_nocuda) (push) Successful in 32m48s
Build Packages / build:rpm (rocky9_sls9) (push) Successful in 35m27s
Build Packages / Generate python client (push) Successful in 25s
Build Packages / build:rpm (rocky9) (push) Successful in 31m59s
Build Packages / Create release (push) Skipped
Build Packages / Build documentation (push) Successful in 1m36s
Build Packages / build:rpm (ubuntu2204) (push) Successful in 24m8s
Build Packages / XDS test (neggia plugin) (push) Successful in 17m46s
Build Packages / build:rpm (ubuntu2404) (push) Successful in 21m36s
Build Packages / XDS test (durin plugin) (push) Successful in 19m40s
Build Packages / XDS test (JFJoch plugin) (push) Successful in 19m38s
Build Packages / DIALS test (push) Successful in 26m30s
1126 lines
51 KiB
C++
1126 lines
51 KiB
C++
// SPDX-FileCopyrightText: 2025 Filip Leonarski, Paul Scherrer Institute <filip.leonarski@psi.ch>
|
|
// SPDX-License-Identifier: GPL-3.0-only
|
|
|
|
#include "PixelRefine.h"
|
|
|
|
#include <Eigen/Dense>
|
|
#include <ceres/ceres.h>
|
|
#include <ceres/rotation.h>
|
|
|
|
#include "../geom_refinement/LatticeReduction.h"
|
|
|
|
namespace {
|
|
|
|
// Per-pixel observation, in *corrected* intensity units (solid-angle and
|
|
// polarization correction already folded in, consistently for signal and
|
|
// background). Geometry-independent quantities are precomputed here so that the
|
|
// Ceres cost functor stays cheap.
|
|
struct PixelObs {
|
|
double x, y; // detector pixel coordinate
|
|
double Iobs; // corrected pixel value (signal + background)
|
|
double Ibkg; // corrected background estimate (azimuthal bin mean)
|
|
double weight; // 1 / sigma_pixel
|
|
double A_recip; // reciprocal-space area subtended by the pixel (Jacobian)
|
|
double angle_rad; // goniometer angle of this observation
|
|
};
|
|
|
|
// One reflection together with the pixels of its shoebox.
|
|
struct ReflGroup {
|
|
int h, k, l;
|
|
double d;
|
|
double Itrue; // reference intensity (held fixed)
|
|
double R_bw_sq; // bandwidth radial-width^2 contribution (0 = monochromatic)
|
|
double predicted_x, predicted_y;
|
|
std::vector<PixelObs> pixels;
|
|
};
|
|
|
|
double SafeInv(double x, double fallback) {
|
|
if (!std::isfinite(x) || std::fabs(x) < 1e-30)
|
|
return fallback;
|
|
return 1.0 / x;
|
|
}
|
|
|
|
// Per-pixel: map a detector pixel through the current geometry into the
|
|
// reference reciprocal frame. Cheap (a few trig + one rotation); depends on the
|
|
// pixel and the detector geometry, not on the lattice.
|
|
template<typename T>
|
|
void ObservedRecip(const T *beam, const T *distance_mm, const T *detector_rot,
|
|
const T *rotation_axis, double obs_x, double obs_y,
|
|
double pixel_size, double inv_lambda, double angle_rad,
|
|
Eigen::Matrix<T, 3, 1> &e_obs_recip) {
|
|
// PyFAI convention (left-handed for rot1/rot2): rot3 = 0 assumed.
|
|
const T c1 = ceres::cos(detector_rot[0]);
|
|
const T s1 = ceres::sin(detector_rot[0]);
|
|
const T c2 = ceres::cos(detector_rot[1]);
|
|
const T s2 = ceres::sin(detector_rot[1]);
|
|
|
|
const T det_x = (T(obs_x) - beam[0]) * T(pixel_size);
|
|
const T det_y = (T(obs_y) - beam[1]) * T(pixel_size);
|
|
const T det_z = T(distance_mm[0]);
|
|
|
|
const T t1_x = c1 * det_x + s1 * det_z;
|
|
const T t1_y = det_y;
|
|
const T t1_z = -s1 * det_x + c1 * det_z;
|
|
|
|
const T x = t1_x;
|
|
const T y = c2 * t1_y + s2 * t1_z;
|
|
const T z = -s2 * t1_y + c2 * t1_z;
|
|
|
|
const T inv_norm = T(1) / ceres::sqrt(x * x + y * y + z * z);
|
|
|
|
T recip_raw[3] = {
|
|
x * inv_norm * T(inv_lambda),
|
|
y * inv_norm * T(inv_lambda),
|
|
(z * inv_norm - T(1.0)) * T(inv_lambda)
|
|
};
|
|
const T aa_back[3] = {
|
|
T(angle_rad) * rotation_axis[0],
|
|
T(angle_rad) * rotation_axis[1],
|
|
T(angle_rad) * rotation_axis[2]
|
|
};
|
|
T recip_obs[3];
|
|
ceres::AngleAxisRotatePoint(aa_back, recip_raw, recip_obs);
|
|
e_obs_recip = Eigen::Matrix<T, 3, 1>(recip_obs[0], recip_obs[1], recip_obs[2]);
|
|
}
|
|
|
|
// Per-reflection: predicted node g_hkl, |g_hkl|^2, and the Ewald-sphere normal.
|
|
// This is the expensive part (symmetry-aware B matrix, three rotations, cross
|
|
// products) - it depends only on the lattice (p0,p1,p2) and hkl, so for a whole
|
|
// shoebox it can be computed once. Convention identical to XtalOptimizer.
|
|
template<typename T>
|
|
bool PredictedNode(const T *p0, const T *p1, const T *p2,
|
|
double exp_h, double exp_k, double exp_l,
|
|
gemmi::CrystalSystem symmetry, double inv_lambda,
|
|
Eigen::Matrix<T, 3, 1> &e_pred_recip,
|
|
Eigen::Matrix<T, 3, 1> &n_radial, T &q_sq) {
|
|
Eigen::Matrix<T, 3, 1> e_uc_len = Eigen::Matrix<T, 3, 1>::Zero();
|
|
Eigen::Matrix<T, 3, 3> Bmat = Eigen::Matrix<T, 3, 3>::Identity();
|
|
|
|
if (symmetry == gemmi::CrystalSystem::Hexagonal) {
|
|
e_uc_len << p1[0], p1[0], p1[2];
|
|
Bmat(0, 1) = T(-0.5);
|
|
Bmat(1, 1) = T(sqrt(3.0) / 2.0);
|
|
} else if (symmetry == gemmi::CrystalSystem::Orthorhombic) {
|
|
e_uc_len << p1[0], p1[1], p1[2];
|
|
} else if (symmetry == gemmi::CrystalSystem::Tetragonal) {
|
|
e_uc_len << p1[0], p1[0], p1[2];
|
|
} else if (symmetry == gemmi::CrystalSystem::Cubic) {
|
|
e_uc_len << p1[0], p1[0], p1[0];
|
|
} else if (symmetry == gemmi::CrystalSystem::Monoclinic) {
|
|
e_uc_len << p1[0], p1[1], p1[2];
|
|
Bmat(0, 2) = ceres::cos(p2[0]);
|
|
Bmat(2, 2) = ceres::sin(p2[0]);
|
|
} else {
|
|
const T ca = ceres::cos(p2[0]);
|
|
const T cb = ceres::cos(p2[1]);
|
|
const T cg = ceres::cos(p2[2]);
|
|
const T sg = ceres::sin(p2[2]);
|
|
|
|
e_uc_len << p1[0], p1[1], p1[2];
|
|
|
|
Bmat(0, 1) = cg;
|
|
Bmat(1, 1) = sg;
|
|
|
|
const T cx = cb;
|
|
const T cy = (ca - cb * cg) / sg;
|
|
const T v = T(1) - cx * cx - cy * cy;
|
|
const T cz = (v >= T(0)) ? ceres::sqrt(v) : T(0);
|
|
|
|
Bmat(0, 2) = cx;
|
|
Bmat(1, 2) = cy;
|
|
Bmat(2, 2) = cz;
|
|
}
|
|
|
|
const T L0 = e_uc_len[0];
|
|
const T L1 = e_uc_len[1];
|
|
const T L2 = e_uc_len[2];
|
|
|
|
T col0_unrot[3] = {Bmat(0, 0) * L0, Bmat(1, 0) * L0, Bmat(2, 0) * L0};
|
|
T col1_unrot[3] = {Bmat(0, 1) * L1, Bmat(1, 1) * L1, Bmat(2, 1) * L1};
|
|
T col2_unrot[3] = {Bmat(0, 2) * L2, Bmat(1, 2) * L2, Bmat(2, 2) * L2};
|
|
|
|
T col0_rot[3], col1_rot[3], col2_rot[3];
|
|
ceres::AngleAxisRotatePoint(p0, col0_unrot, col0_rot);
|
|
ceres::AngleAxisRotatePoint(p0, col1_unrot, col1_rot);
|
|
ceres::AngleAxisRotatePoint(p0, col2_unrot, col2_rot);
|
|
|
|
const Eigen::Matrix<T, 3, 1> A(col0_rot[0], col0_rot[1], col0_rot[2]);
|
|
const Eigen::Matrix<T, 3, 1> Bv(col1_rot[0], col1_rot[1], col1_rot[2]);
|
|
const Eigen::Matrix<T, 3, 1> C(col2_rot[0], col2_rot[1], col2_rot[2]);
|
|
|
|
const Eigen::Matrix<T, 3, 1> BxC = Bv.cross(C);
|
|
const Eigen::Matrix<T, 3, 1> CxA = C.cross(A);
|
|
const Eigen::Matrix<T, 3, 1> AxB = A.cross(Bv);
|
|
|
|
const T Vol = A.dot(BxC);
|
|
if (ceres::abs(Vol) < T(1e-12))
|
|
return false;
|
|
const T invV = T(1) / Vol;
|
|
|
|
e_pred_recip = (BxC * T(exp_h) + CxA * T(exp_k) + AxB * T(exp_l)) * invV;
|
|
q_sq = e_pred_recip.squaredNorm();
|
|
|
|
// Ewald sphere centre at -k_i = (0,0,-inv_lambda); radial normal at g_hkl.
|
|
const Eigen::Matrix<T, 3, 1> S_pred(
|
|
e_pred_recip[0],
|
|
e_pred_recip[1],
|
|
e_pred_recip[2] + T(inv_lambda));
|
|
const T S_pred_norm = S_pred.norm();
|
|
if (S_pred_norm < T(1e-10))
|
|
return false;
|
|
|
|
n_radial = S_pred / S_pred_norm;
|
|
return true;
|
|
}
|
|
|
|
} // namespace
|
|
|
|
// ---------------------------------------------------------------------------
|
|
// Cost functor
|
|
//
|
|
// I_pred(pixel) = G * Itrue * B_term * P_radial * P_tangential + I_bkg
|
|
//
|
|
// B_term = exp(-B |q|^2 / 4) (Debye-Waller)
|
|
// P_radial = exp(-eps_r^2 / R0_eff^2) (partiality: fraction of
|
|
// the mosaic blob on the
|
|
// Ewald sphere; <= 1)
|
|
// P_tangential = A_recip/(pi R1^2) * exp(-eps_t^2/R1^2)(spatial profile on the
|
|
// detector, normalized so
|
|
// that sum over pixels ~ 1)
|
|
//
|
|
// The tangential factor is what makes this "profile fitting": summing
|
|
// I_pred - I_bkg over the shoebox reproduces G * Itrue * B_term * P_radial.
|
|
// The 1/(pi R1^2) normalization is the missing piece that decouples the profile
|
|
// width R1 from the overall scale G.
|
|
//
|
|
// X-ray bandwidth: a spread in lambda is a spread in the Ewald-sphere radius,
|
|
// i.e. a purely *radial* thickening of the shell. It adds (in quadrature) a
|
|
// resolution-dependent term to the radial width:
|
|
// R0_eff^2 = R0^2 + R_bw^2 , R_bw^2 = (b*lambda)^2 / (2 d^4)
|
|
// where b = relative bandwidth (sigma of dlambda/lambda). R_bw grows like 1/d^2,
|
|
// so bandwidth leaves low-resolution spots sharp and smears high-resolution ones
|
|
// radially - the pink-beam/DMM signature. R_bw_sq is a fixed per-reflection
|
|
// constant (b is known), so R0 keeps meaning "intrinsic" width (mosaic +
|
|
// divergence + beam). b = 0 makes R_bw = 0: a monochromatic no-op.
|
|
// ---------------------------------------------------------------------------
|
|
struct PixelResidual {
|
|
PixelResidual(const PixelObs &obs, double Itrue,
|
|
double lambda, double pixel_size,
|
|
double exp_h, double exp_k, double exp_l,
|
|
double R_bw_sq,
|
|
gemmi::CrystalSystem symmetry)
|
|
: Itrue(Itrue), Iobs(obs.Iobs), Ibkg(obs.Ibkg), weight(obs.weight),
|
|
A_recip(obs.A_recip), obs_x(obs.x), obs_y(obs.y),
|
|
inv_lambda(1.0 / lambda), pixel_size(pixel_size),
|
|
exp_h(exp_h), exp_k(exp_k), exp_l(exp_l),
|
|
R_bw_sq(R_bw_sq), angle_rad(obs.angle_rad), symmetry(symmetry) {
|
|
if (std::fabs(lambda) < 1e-6)
|
|
throw JFJochException(JFJochExceptionCategory::InputParameterInvalid,
|
|
"Lambda cannot be close to zero");
|
|
}
|
|
|
|
// Maps a detector pixel through the current geometry + lattice into the
|
|
// reference reciprocal frame and returns:
|
|
// q_sq = |g_hkl|^2 (predicted node, for B-factor)
|
|
// eps_radial = deviation along Ewald normal (partiality direction)
|
|
// eps_tang_sq = squared deviation in the detector-tangential plane (profile)
|
|
template<typename T>
|
|
bool GeometryTerms(const T *const beam,
|
|
const T *const distance_mm,
|
|
const T *const detector_rot,
|
|
const T *const rotation_axis,
|
|
const T *const p0,
|
|
const T *const p1,
|
|
const T *const p2,
|
|
T &q_sq, T &eps_radial, T &eps_tang_sq) const {
|
|
Eigen::Matrix<T, 3, 1> e_obs_recip;
|
|
ObservedRecip(beam, distance_mm, detector_rot, rotation_axis,
|
|
obs_x, obs_y, pixel_size, inv_lambda, angle_rad, e_obs_recip);
|
|
|
|
Eigen::Matrix<T, 3, 1> e_pred_recip, n_radial;
|
|
if (!PredictedNode(p0, p1, p2, exp_h, exp_k, exp_l, symmetry, inv_lambda,
|
|
e_pred_recip, n_radial, q_sq))
|
|
return false;
|
|
|
|
const Eigen::Matrix<T, 3, 1> delta_q = e_obs_recip - e_pred_recip;
|
|
eps_radial = delta_q.dot(n_radial);
|
|
eps_tang_sq = (delta_q - eps_radial * n_radial).squaredNorm();
|
|
return true;
|
|
}
|
|
|
|
// Assembles the full model intensity for the pixel from the geometry terms.
|
|
template<typename T>
|
|
bool Model(const T *const beam, const T *const distance_mm,
|
|
const T *const detector_rot, const T *const rotation_axis,
|
|
const T *const p0, const T *const p1, const T *const p2,
|
|
const T *const scale_factor, const T *const B, const T *const R,
|
|
T &Ipred) const {
|
|
T q_sq, eps_radial, eps_tang_sq;
|
|
if (!GeometryTerms(beam, distance_mm, detector_rot, rotation_axis,
|
|
p0, p1, p2, q_sq, eps_radial, eps_tang_sq))
|
|
return false;
|
|
|
|
if (R[0] < T(1e-10) || R[1] < T(1e-10))
|
|
return false;
|
|
|
|
const T B_term = ceres::exp(-B[0] * q_sq / T(4.0));
|
|
|
|
// Separable Gaussian spot model:
|
|
// radial P_r(e) = exp(-e^2/R0_eff^2) (peak-normalized, in (0,1])
|
|
// tangent g_t(e) = exp(-|e|^2/R1^2) / (pi R1^2) [1/A^-2]
|
|
// The pixel captures the fraction g_t * A_recip of the tangential profile
|
|
// (A_recip = reciprocal area the pixel subtends; sum over shoebox ~ 1).
|
|
// The radial factor is the still-image partiality (how far the reflection
|
|
// sits from the Ewald sphere); the overall scale is carried by the free G.
|
|
//
|
|
// IMPORTANT: the radial factor MUST use the same convention here as the
|
|
// extraction's `partiality` (peak-normalized), otherwise image_scale_corr
|
|
// = 1/(partiality*G*B) does not invert the model and a leftover, R0_eff-
|
|
// dependent (hence resolution-dependent) factor biases the intensities.
|
|
// R0_eff folds in the energy-bandwidth broadening via R_bw_sq.
|
|
const T R0_eff_sq = R[0] * R[0] + T(R_bw_sq);
|
|
const T P_radial = ceres::exp(-eps_radial * eps_radial / R0_eff_sq);
|
|
const T P_tang = T(A_recip) * ceres::exp(-eps_tang_sq / (R[1] * R[1]))
|
|
/ (T(M_PI) * R[1] * R[1]);
|
|
|
|
const T signal = scale_factor[0] * T(Itrue) * B_term * P_radial * P_tang;
|
|
Ipred = signal + T(Ibkg);
|
|
return true;
|
|
}
|
|
|
|
template<typename T>
|
|
bool operator()(const T *const beam,
|
|
const T *const distance_mm,
|
|
const T *const detector_rot,
|
|
const T *const rotation_axis,
|
|
const T *const p0,
|
|
const T *const p1,
|
|
const T *const p2,
|
|
const T *const scale_factor,
|
|
const T *const B,
|
|
const T *const R,
|
|
T *residual) const {
|
|
T Ipred;
|
|
if (!Model(beam, distance_mm, detector_rot, rotation_axis,
|
|
p0, p1, p2, scale_factor, B, R, Ipred))
|
|
return false;
|
|
|
|
residual[0] = (Ipred - T(Iobs)) * T(weight);
|
|
return true;
|
|
}
|
|
|
|
const double Itrue, Iobs, Ibkg, weight, A_recip;
|
|
const double obs_x, obs_y;
|
|
const double inv_lambda;
|
|
const double pixel_size;
|
|
const double exp_h, exp_k, exp_l;
|
|
const double R_bw_sq; // bandwidth radial-width^2 contribution (0 = monochromatic)
|
|
const double angle_rad;
|
|
gemmi::CrystalSystem symmetry;
|
|
};
|
|
|
|
// ---------------------------------------------------------------------------
|
|
// Per-shoebox cost functor
|
|
//
|
|
// One residual block per reflection emitting N residuals (one per shoebox pixel).
|
|
// The expensive per-reflection geometry (PredictedNode: symmetry-aware B matrix,
|
|
// three rotations, cross products) is computed ONCE; only the cheap per-pixel
|
|
// ObservedRecip + Gaussian profile run in the pixel loop. This is identical in
|
|
// value to the old one-block-per-pixel formulation but ~(pixels-per-shoebox)x
|
|
// fewer evaluations of the costly node computation. Uses the same shared helpers
|
|
// (and hence the same conventions) as PixelResidual.
|
|
// ---------------------------------------------------------------------------
|
|
struct ShoeboxResidual {
|
|
ShoeboxResidual(const ReflGroup &g, double lambda, double pixel_size,
|
|
gemmi::CrystalSystem symmetry)
|
|
: pixels(g.pixels), Itrue(g.Itrue), R_bw_sq(g.R_bw_sq),
|
|
exp_h(g.h), exp_k(g.k), exp_l(g.l),
|
|
inv_lambda(1.0 / lambda), pixel_size(pixel_size),
|
|
angle_rad(g.pixels.empty() ? 0.0 : g.pixels.front().angle_rad),
|
|
symmetry(symmetry) {}
|
|
|
|
template<typename T>
|
|
bool operator()(const T *const *params, T *residual) const {
|
|
// Parameter blocks (order matches AddParameterBlock in Run):
|
|
// 0 beam[2] 1 distance[1] 2 detector_rot[2] 3 rotation_axis[3]
|
|
// 4 p0[3] 5 p1[3] 6 p2[3] 7 scale[1] 8 B[1] 9 R[2]
|
|
const T *beam = params[0];
|
|
const T *distance_mm = params[1];
|
|
const T *detector_rot = params[2];
|
|
const T *rotation_axis = params[3];
|
|
const T *p0 = params[4];
|
|
const T *p1 = params[5];
|
|
const T *p2 = params[6];
|
|
const T *scale_factor = params[7];
|
|
const T *B = params[8];
|
|
const T *R = params[9];
|
|
|
|
if (R[0] < T(1e-10) || R[1] < T(1e-10))
|
|
return false;
|
|
|
|
// --- per-reflection: computed once ---------------------------------
|
|
Eigen::Matrix<T, 3, 1> e_pred_recip, n_radial;
|
|
T q_sq;
|
|
if (!PredictedNode(p0, p1, p2, exp_h, exp_k, exp_l, symmetry, inv_lambda,
|
|
e_pred_recip, n_radial, q_sq))
|
|
return false;
|
|
|
|
const T B_term = ceres::exp(-B[0] * q_sq / T(4.0));
|
|
const T R0_eff_sq = R[0] * R[0] + T(R_bw_sq);
|
|
|
|
// --- per-pixel loop -------------------------------------------------
|
|
for (size_t i = 0; i < pixels.size(); ++i) {
|
|
const PixelObs &obs = pixels[i];
|
|
|
|
Eigen::Matrix<T, 3, 1> e_obs_recip;
|
|
ObservedRecip(beam, distance_mm, detector_rot, rotation_axis,
|
|
obs.x, obs.y, pixel_size, inv_lambda, angle_rad, e_obs_recip);
|
|
|
|
const Eigen::Matrix<T, 3, 1> delta_q = e_obs_recip - e_pred_recip;
|
|
const T eps_radial = delta_q.dot(n_radial);
|
|
const T eps_tang_sq = (delta_q - eps_radial * n_radial).squaredNorm();
|
|
|
|
const T P_radial = ceres::exp(-eps_radial * eps_radial / R0_eff_sq);
|
|
const T P_tang = T(obs.A_recip) * ceres::exp(-eps_tang_sq / (R[1] * R[1]))
|
|
/ (T(M_PI) * R[1] * R[1]);
|
|
|
|
const T signal = scale_factor[0] * T(Itrue) * B_term * P_radial * P_tang;
|
|
const T Ipred = signal + T(obs.Ibkg);
|
|
residual[i] = (Ipred - T(obs.Iobs)) * T(obs.weight);
|
|
}
|
|
return true;
|
|
}
|
|
|
|
std::vector<PixelObs> pixels;
|
|
const double Itrue, R_bw_sq;
|
|
const double exp_h, exp_k, exp_l;
|
|
const double inv_lambda, pixel_size, angle_rad;
|
|
gemmi::CrystalSystem symmetry;
|
|
};
|
|
|
|
PixelRefine::PixelRefine(const DiffractionExperiment &experiment,
|
|
const AzimuthalIntegrationMapping &mapping,
|
|
const std::vector<MergedReflection> &reference)
|
|
: mapping(mapping),
|
|
xpixel(experiment.GetXPixelsNum()),
|
|
ypixel(experiment.GetYPixelsNum()),
|
|
experiment(experiment),
|
|
hkl_key_generator(experiment.GetScalingSettings().GetMergeFriedel(),
|
|
experiment.GetSpaceGroupNumber().value_or(1)) {
|
|
for (const auto &ref: reference)
|
|
reference_data[hkl_key_generator(ref)] = ref.I;
|
|
}
|
|
|
|
void PixelRefine::BuildParameterBlocks(const PixelRefineData &data,
|
|
double beam[2], double &dist_mm,
|
|
double detector_rot[2], double rot_vec[3],
|
|
double latt_vec0[3], double latt_vec1[3], double latt_vec2[3]) const {
|
|
beam[0] = data.geom.GetBeamX_pxl();
|
|
beam[1] = data.geom.GetBeamY_pxl();
|
|
dist_mm = data.geom.GetDetectorDistance_mm();
|
|
detector_rot[0] = data.geom.GetPoniRot1_rad();
|
|
detector_rot[1] = data.geom.GetPoniRot2_rad();
|
|
rot_vec[0] = 1.0; rot_vec[1] = 0.0; rot_vec[2] = 0.0;
|
|
if (auto axis = data.geom.GetRotation()) {
|
|
rot_vec[0] = axis->GetAxis().x;
|
|
rot_vec[1] = axis->GetAxis().y;
|
|
rot_vec[2] = axis->GetAxis().z;
|
|
}
|
|
|
|
for (int i = 0; i < 3; ++i)
|
|
latt_vec0[i] = latt_vec1[i] = latt_vec2[i] = 0.0;
|
|
|
|
double beta = data.latt.GetUnitCell().beta;
|
|
switch (data.crystal_system) {
|
|
case gemmi::CrystalSystem::Orthorhombic:
|
|
LatticeToRodriguesAndLengths_GS(data.latt, latt_vec0, latt_vec1);
|
|
break;
|
|
case gemmi::CrystalSystem::Tetragonal:
|
|
LatticeToRodriguesAndLengths_GS(data.latt, latt_vec0, latt_vec1);
|
|
latt_vec1[0] = (latt_vec1[0] + latt_vec1[1]) / 2.0;
|
|
break;
|
|
case gemmi::CrystalSystem::Cubic:
|
|
LatticeToRodriguesAndLengths_GS(data.latt, latt_vec0, latt_vec1);
|
|
latt_vec1[0] = (latt_vec1[0] + latt_vec1[1] + latt_vec1[2]) / 3.0;
|
|
break;
|
|
case gemmi::CrystalSystem::Hexagonal:
|
|
LatticeToRodriguesAndLengths_Hex(data.latt, latt_vec0, latt_vec1);
|
|
break;
|
|
case gemmi::CrystalSystem::Monoclinic:
|
|
LatticeToRodriguesLengthsBeta_Mono(data.latt, latt_vec0, latt_vec1, beta);
|
|
latt_vec2[0] = beta;
|
|
break;
|
|
default: {
|
|
LatticeToRodriguesAndLengths_GS(data.latt, latt_vec0, latt_vec1);
|
|
const auto uc = data.latt.GetUnitCell();
|
|
latt_vec2[0] = uc.alpha * M_PI / 180.0;
|
|
latt_vec2[1] = uc.beta * M_PI / 180.0;
|
|
latt_vec2[2] = uc.gamma * M_PI / 180.0;
|
|
break;
|
|
}
|
|
}
|
|
}
|
|
|
|
template<class T>
|
|
void PixelRefine::Run(const T *image,
|
|
const AzimuthalIntegrationProfile &profile,
|
|
BraggPrediction &prediction,
|
|
PixelRefineData &data) {
|
|
data.solved = false;
|
|
data.reflections.clear();
|
|
|
|
const double lambda = data.geom.GetWavelength_A();
|
|
const double pixel_size = data.geom.GetPixelSize_mm();
|
|
|
|
const BraggPredictionSettings settings_prediction{
|
|
.high_res_A = experiment.GetBraggIntegrationSettings().GetDMinLimit_A(),
|
|
.max_hkl = 100,
|
|
.centering = data.centering,
|
|
.bandwidth_sigma = static_cast<float>(data.bandwidth) // relative Δλ/λ sigma
|
|
};
|
|
|
|
const auto azim_result = profile.GetResult();
|
|
const auto azim_std = profile.GetStd();
|
|
const auto &pixel_to_bin = mapping.GetPixelToBin();
|
|
const auto &corrections = mapping.Corrections();
|
|
// pixel_to_bin stores the *full* bin index (azimuthal_sector * q_bins + q_bin),
|
|
// so the valid range is the total number of bins, i.e. the profile size - NOT
|
|
// GetAzimuthalBinCount() (which is only the number of azimuthal sectors).
|
|
const int total_bin_count = static_cast<int>(azim_result.size());
|
|
|
|
const double angle_rad = data.angle_deg * M_PI / 180.0;
|
|
const int radius = data.shoebox_radius;
|
|
|
|
// Exact reciprocal-space area a 1x1 pixel subtends, |dq/dx x dq/dy|, via
|
|
// finite differences of the detector->reciprocal map. This is the Jacobian
|
|
// between the curved Ewald-sphere sampling and flat reciprocal space, and it
|
|
// is exactly the geometric factor that plays the role of the Lorentz factor
|
|
// for stills: where the sphere grazes reciprocal space obliquely, a pixel
|
|
// covers more reciprocal volume and the captured fraction grows. It tracks
|
|
// the refined geometry because it reads the current data.geom each iteration.
|
|
auto recip_area = [&](double x, double y) -> double {
|
|
const Coord qx = data.geom.DetectorToRecip(x + 0.5, y) - data.geom.DetectorToRecip(x - 0.5, y);
|
|
const Coord qy = data.geom.DetectorToRecip(x, y + 0.5) - data.geom.DetectorToRecip(x, y - 0.5);
|
|
return (qx % qy).Length();
|
|
};
|
|
|
|
// Bandwidth radial-width^2 (in the code's R = sqrt(2)*sigma convention):
|
|
// R_bw^2 = (b*lambda)^2 / (2 d^4), b = relative bandwidth (sigma).
|
|
// A fixed per-reflection constant; data.bandwidth == 0 -> monochromatic no-op.
|
|
const double bw = data.bandwidth;
|
|
auto bandwidth_radial_sq = [&](double d) -> double {
|
|
if (bw <= 0.0 || d <= 0.0)
|
|
return 0.0;
|
|
const double bl = bw * lambda;
|
|
return bl * bl / (2.0 * d * d * d * d);
|
|
};
|
|
|
|
// Mutable experiment whose geometry is re-synced from the refined data.geom
|
|
// before each prediction, so shoeboxes track the refined geometry/cell.
|
|
DiffractionExperiment exp_iter = experiment;
|
|
|
|
// State retained after the loop for the final reflection extraction.
|
|
std::vector<ReflGroup> groups;
|
|
double beam[2] = {0, 0};
|
|
double dist_mm = data.geom.GetDetectorDistance_mm();
|
|
double detector_rot[2] = {0, 0};
|
|
double rot_vec[3] = {1.0, 0.0, 0.0};
|
|
double latt_vec0[3] = {0, 0, 0}; // orientation (Rodrigues)
|
|
double latt_vec1[3] = {0, 0, 0}; // lengths
|
|
double latt_vec2[3] = {0, 0, 0}; // angles (rad)
|
|
|
|
const bool eval_only = (data.max_iterations <= 0);
|
|
const int n_iter = std::max(1, data.max_iterations);
|
|
for (int iter = 0; iter < n_iter; ++iter) {
|
|
// ---- 1. Re-sync prediction geometry from the (refined) model ----------
|
|
exp_iter.BeamX_pxl(data.geom.GetBeamX_pxl())
|
|
.BeamY_pxl(data.geom.GetBeamY_pxl())
|
|
.DetectorDistance_mm(data.geom.GetDetectorDistance_mm())
|
|
.PoniRot1_rad(data.geom.GetPoniRot1_rad())
|
|
.PoniRot2_rad(data.geom.GetPoniRot2_rad());
|
|
|
|
const int nrefl = prediction.Calc(exp_iter, data.latt, settings_prediction);
|
|
|
|
// ---- 2. Collect per-reflection shoebox pixels -------------------------
|
|
// GetReflections() returns the full pre-sized buffer; only the first
|
|
// nrefl entries are valid for this image (the rest are stale/zeroed).
|
|
groups.clear();
|
|
const auto &predicted = prediction.GetReflections();
|
|
for (int ri = 0; ri < nrefl; ++ri) {
|
|
const auto &refl = predicted[ri];
|
|
const auto hkl = hkl_key_generator(refl);
|
|
if (!reference_data.contains(hkl))
|
|
continue;
|
|
|
|
ReflGroup g;
|
|
g.h = refl.h;
|
|
g.k = refl.k;
|
|
g.l = refl.l;
|
|
g.d = refl.d;
|
|
g.Itrue = reference_data[hkl];
|
|
g.R_bw_sq = bandwidth_radial_sq(refl.d);
|
|
g.predicted_x = refl.predicted_x;
|
|
g.predicted_y = refl.predicted_y;
|
|
|
|
const int min_y = std::max<int>(refl.predicted_y - radius, 0);
|
|
const int max_y = std::min<int>(refl.predicted_y + radius, ypixel - 1);
|
|
const int min_x = std::max<int>(refl.predicted_x - radius, 0);
|
|
const int max_x = std::min<int>(refl.predicted_x + radius, xpixel - 1);
|
|
|
|
for (int y = min_y; y <= max_y; ++y) {
|
|
for (int x = min_x; x <= max_x; ++x) {
|
|
const size_t npixel = xpixel * y + x;
|
|
const int azim_bin = pixel_to_bin[npixel];
|
|
|
|
// Skip pixels not mapped to a bin or carrying a sentinel
|
|
// (masked / saturated) value. We assume the pixel mask is
|
|
// already applied upstream.
|
|
if (azim_bin >= total_bin_count)
|
|
continue;
|
|
if (image[npixel] == std::numeric_limits<T>::max())
|
|
continue;
|
|
if (std::is_signed_v<T> && (image[npixel] == std::numeric_limits<T>::min()))
|
|
continue;
|
|
|
|
const double correction = corrections[npixel];
|
|
const double Ibkg = azim_result[azim_bin]; // already in corrected units
|
|
const double Ibkg_sigma = azim_std[azim_bin];
|
|
const double raw = static_cast<double>(image[npixel]);
|
|
const double Iobs = raw * correction;
|
|
|
|
// Per-pixel variance: Poisson noise of the corrected counts
|
|
// (var(c*N) = c^2 * N = c * Iobs) plus the background spread.
|
|
double var = correction * std::max(Iobs, 0.0) + Ibkg_sigma * Ibkg_sigma;
|
|
if (!(var > 1.0))
|
|
var = 1.0;
|
|
|
|
PixelObs obs{
|
|
.x = static_cast<double>(x),
|
|
.y = static_cast<double>(y),
|
|
.Iobs = Iobs,
|
|
.Ibkg = Ibkg,
|
|
.weight = 1.0 / std::sqrt(var),
|
|
.A_recip = recip_area(x, y),
|
|
.angle_rad = angle_rad
|
|
};
|
|
g.pixels.push_back(obs);
|
|
}
|
|
}
|
|
|
|
if (!g.pixels.empty())
|
|
groups.push_back(std::move(g));
|
|
}
|
|
|
|
if (groups.empty())
|
|
return;
|
|
|
|
// ---- 3. Set up parameter blocks (geometry part mirrors XtalOptimizer) -
|
|
BuildParameterBlocks(data, beam, dist_mm, detector_rot, rot_vec,
|
|
latt_vec0, latt_vec1, latt_vec2);
|
|
|
|
// ---- 4. Build the problem ---------------------------------------------
|
|
// One residual block per shoebox (N residuals), so the expensive
|
|
// per-reflection node geometry is evaluated once per reflection instead
|
|
// of once per pixel.
|
|
ceres::Problem problem;
|
|
size_t residual_pixels = 0;
|
|
for (const auto &g : groups) {
|
|
auto *cost = new ceres::DynamicAutoDiffCostFunction<ShoeboxResidual>(
|
|
new ShoeboxResidual(g, lambda, pixel_size, data.crystal_system));
|
|
cost->AddParameterBlock(2); // beam
|
|
cost->AddParameterBlock(1); // distance
|
|
cost->AddParameterBlock(2); // detector_rot
|
|
cost->AddParameterBlock(3); // rotation_axis
|
|
cost->AddParameterBlock(3); // p0 (orientation)
|
|
cost->AddParameterBlock(3); // p1 (lengths)
|
|
cost->AddParameterBlock(3); // p2 (angles)
|
|
cost->AddParameterBlock(1); // scale G
|
|
cost->AddParameterBlock(1); // B
|
|
cost->AddParameterBlock(2); // R
|
|
cost->SetNumResiduals(static_cast<int>(g.pixels.size()));
|
|
// No robust loss here: a per-block (whole-shoebox) Huber would act on
|
|
// the sum of ~N squared residuals and mis-scale, unlike the previous
|
|
// per-pixel Huber. Per-pixel sigma weighting is retained; per-pixel
|
|
// outlier rejection (zingers) is a TODO if needed.
|
|
problem.AddResidualBlock(cost, nullptr,
|
|
beam, &dist_mm, detector_rot, rot_vec,
|
|
latt_vec0, latt_vec1, latt_vec2,
|
|
&data.scale_factor, &data.B_factor, data.R);
|
|
residual_pixels += g.pixels.size();
|
|
}
|
|
data.residual_count = residual_pixels;
|
|
|
|
// ---- 5. Constrain / bound parameter blocks ----------------------------
|
|
if (!data.refine_orientation)
|
|
problem.SetParameterBlockConstant(latt_vec0);
|
|
|
|
if (!data.refine_unit_cell) {
|
|
problem.SetParameterBlockConstant(latt_vec1);
|
|
problem.SetParameterBlockConstant(latt_vec2);
|
|
} else {
|
|
for (int i = 0; i < 3; ++i) {
|
|
problem.SetParameterLowerBound(latt_vec1, i, 5.0);
|
|
problem.SetParameterUpperBound(latt_vec1, i, 1000.0);
|
|
}
|
|
if (data.crystal_system != gemmi::CrystalSystem::Monoclinic &&
|
|
data.crystal_system != gemmi::CrystalSystem::Triclinic)
|
|
problem.SetParameterBlockConstant(latt_vec2);
|
|
}
|
|
|
|
if (!data.refine_beam_center)
|
|
problem.SetParameterBlockConstant(beam);
|
|
|
|
if (!data.refine_distance) {
|
|
problem.SetParameterBlockConstant(&dist_mm);
|
|
} else {
|
|
problem.SetParameterLowerBound(&dist_mm, 0, dist_mm * 0.9);
|
|
problem.SetParameterUpperBound(&dist_mm, 0, dist_mm * 1.1);
|
|
}
|
|
|
|
if (!data.refine_detector_angles) {
|
|
problem.SetParameterBlockConstant(detector_rot);
|
|
} else {
|
|
const double rng = 3.0 / 180.0 * M_PI;
|
|
for (int i = 0; i < 2; ++i) {
|
|
problem.SetParameterLowerBound(detector_rot, i, detector_rot[i] - rng);
|
|
problem.SetParameterUpperBound(detector_rot, i, detector_rot[i] + rng);
|
|
}
|
|
}
|
|
|
|
if (!data.refine_rotation_axis)
|
|
problem.SetParameterBlockConstant(rot_vec);
|
|
|
|
if (data.refine_scale)
|
|
problem.SetParameterLowerBound(&data.scale_factor, 0, 0.0);
|
|
else
|
|
problem.SetParameterBlockConstant(&data.scale_factor);
|
|
|
|
if (!data.refine_B)
|
|
problem.SetParameterBlockConstant(&data.B_factor);
|
|
|
|
if (data.refine_R) {
|
|
problem.SetParameterLowerBound(data.R, 0, 1e-5);
|
|
problem.SetParameterLowerBound(data.R, 1, 1e-5);
|
|
} else {
|
|
problem.SetParameterBlockConstant(data.R);
|
|
}
|
|
|
|
// ---- 6. Solve (or, for max_iterations<=0, just evaluate the cost) -----
|
|
// Evaluate-only is the live-residual path: it reports the current cost
|
|
// without moving any parameter, so a UI can show how good the present
|
|
// R0/R1/bandwidth/geometry are as the user drags sliders.
|
|
if (eval_only) {
|
|
double cost = 0.0;
|
|
problem.Evaluate(ceres::Problem::EvaluateOptions(), &cost, nullptr, nullptr, nullptr);
|
|
data.final_cost = cost;
|
|
data.solved = true;
|
|
} else {
|
|
ceres::Solver::Options options;
|
|
options.linear_solver_type = ceres::DENSE_QR;
|
|
options.minimizer_progress_to_stdout = false;
|
|
options.logging_type = ceres::LoggingType::SILENT;
|
|
options.max_solver_time_in_seconds = data.max_time_s;
|
|
options.num_threads = 1;
|
|
|
|
ceres::Solver::Summary summary;
|
|
ceres::Solve(options, &problem, &summary);
|
|
|
|
data.final_cost = summary.final_cost;
|
|
data.solved = summary.IsSolutionUsable();
|
|
}
|
|
|
|
// ---- 7. Write refined geometry + lattice back into data ---------------
|
|
if (data.refine_beam_center)
|
|
data.geom.BeamX_pxl(beam[0]).BeamY_pxl(beam[1]);
|
|
if (data.refine_distance)
|
|
data.geom.DetectorDistance_mm(dist_mm);
|
|
if (data.refine_detector_angles)
|
|
data.geom.PoniRot1_rad(detector_rot[0]).PoniRot2_rad(detector_rot[1]);
|
|
|
|
if (data.refine_orientation || data.refine_unit_cell) {
|
|
switch (data.crystal_system) {
|
|
case gemmi::CrystalSystem::Orthorhombic:
|
|
data.latt = AngleAxisAndCellToLattice(latt_vec0, latt_vec1, M_PI/2, M_PI/2, M_PI/2);
|
|
break;
|
|
case gemmi::CrystalSystem::Tetragonal:
|
|
latt_vec1[1] = latt_vec1[0];
|
|
data.latt = AngleAxisAndCellToLattice(latt_vec0, latt_vec1, M_PI/2, M_PI/2, M_PI/2);
|
|
break;
|
|
case gemmi::CrystalSystem::Cubic:
|
|
latt_vec1[1] = latt_vec1[0];
|
|
latt_vec1[2] = latt_vec1[0];
|
|
data.latt = AngleAxisAndCellToLattice(latt_vec0, latt_vec1, M_PI/2, M_PI/2, M_PI/2);
|
|
break;
|
|
case gemmi::CrystalSystem::Hexagonal:
|
|
latt_vec1[1] = latt_vec1[0];
|
|
data.latt = AngleAxisAndCellToLattice(latt_vec0, latt_vec1, M_PI/2, M_PI/2, 2.0*M_PI/3.0);
|
|
break;
|
|
case gemmi::CrystalSystem::Monoclinic:
|
|
data.latt = AngleAxisAndCellToLattice(latt_vec0, latt_vec1, M_PI/2, latt_vec2[0], M_PI/2);
|
|
break;
|
|
default:
|
|
data.latt = AngleAxisAndCellToLattice(latt_vec0, latt_vec1,
|
|
latt_vec2[0], latt_vec2[1], latt_vec2[2]);
|
|
break;
|
|
}
|
|
}
|
|
} // predict<->refine iterations
|
|
|
|
// ---- Extract integrated reflections ---------------------------------------
|
|
// Profile fitting gives the recorded amplitude (against the normalized
|
|
// tangential profile P_t):
|
|
// J = sum_p[ P_t,p (Iobs_p - Ibkg_p)/v_p ] / sum_p[ P_t,p^2 / v_p ]
|
|
// ~ G * Itrue * B_term * partiality (recorded intensity)
|
|
// var(J) = 1 / sum_p[ P_t,p^2 / v_p ]
|
|
//
|
|
// Output split (Merge multiplies r.I * image_scale_corr and weights by
|
|
// 1/(sigma*image_scale_corr)^2 - see Merge.cpp):
|
|
// r.I = J / (B_term * partiality) = G * Itrue (B/partiality corrected)
|
|
// r.sigma = sqrt(var(J)) / (B_term * partiality)
|
|
// r.partiality = profile-weighted peak radial factor in (0,1] (Merge filter only)
|
|
// r.image_scale_corr = 1/G (per-image scale ONLY)
|
|
// so r.I * image_scale_corr = Itrue. B and partiality live on the intensity,
|
|
// G lives on image_scale_corr - one clean meaning per field.
|
|
data.reflections.reserve(groups.size());
|
|
for (const auto &g : groups) {
|
|
double num = 0.0, den = 0.0, bkg_sum = 0.0;
|
|
double radial_sum = 0.0, radial_w = 0.0;
|
|
size_t n = 0;
|
|
|
|
for (const auto &obs : g.pixels) {
|
|
PixelResidual pr(obs, 1.0, lambda, pixel_size, g.h, g.k, g.l, g.R_bw_sq, data.crystal_system);
|
|
double q_sq, eps_r, eps_t_sq;
|
|
if (!pr.GeometryTerms(beam, &dist_mm, detector_rot, rot_vec,
|
|
latt_vec0, latt_vec1, latt_vec2, q_sq, eps_r, eps_t_sq))
|
|
continue;
|
|
if (!(data.R[0] > 0.0) || !(data.R[1] > 0.0))
|
|
continue;
|
|
|
|
// Normalized tangential profile (sum over shoebox ~ 1) -> fit weight.
|
|
const double P_t = obs.A_recip * std::exp(-eps_t_sq / (data.R[1] * data.R[1]))
|
|
/ (M_PI * data.R[1] * data.R[1]);
|
|
// Peak-normalized radial factor (the partiality), in (0,1].
|
|
// Bandwidth-broadened radial width, matching the model in Model().
|
|
const double R0_eff_sq = data.R[0] * data.R[0] + g.R_bw_sq;
|
|
const double P_radial = std::exp(-eps_r * eps_r / R0_eff_sq);
|
|
|
|
const double v = SafeInv(obs.weight * obs.weight, 1.0); // pixel variance
|
|
const double signal = obs.Iobs - obs.Ibkg;
|
|
|
|
num += P_t * signal / v;
|
|
den += P_t * P_t / v;
|
|
radial_sum += P_radial * P_t; // weight partiality by the spot core
|
|
radial_w += P_t;
|
|
bkg_sum += obs.Ibkg;
|
|
++n;
|
|
}
|
|
|
|
Reflection r{};
|
|
r.h = g.h;
|
|
r.k = g.k;
|
|
r.l = g.l;
|
|
r.d = static_cast<float>(g.d);
|
|
r.predicted_x = static_cast<float>(g.predicted_x);
|
|
r.predicted_y = static_cast<float>(g.predicted_y);
|
|
r.observed_x = NAN;
|
|
r.observed_y = NAN;
|
|
r.rlp = 1.0f;
|
|
r.partiality = (radial_w > 0.0) ? static_cast<float>(radial_sum / radial_w) : 1.0f;
|
|
|
|
if (den > 0.0 && n > 0) {
|
|
const double I_amp = num / den; // ~ G*Itrue*B_term*partiality
|
|
const double sigma_amp = std::sqrt(1.0 / den);
|
|
const double B_term = std::exp(-data.B_factor / (4.0 * g.d * g.d));
|
|
const double corr = static_cast<double>(r.partiality) * B_term; // B & partiality
|
|
r.bkg = static_cast<float>(bkg_sum / static_cast<double>(n));
|
|
r.observed = true;
|
|
|
|
if (corr > 0.0 && data.scale_factor > 0.0) {
|
|
r.I = static_cast<float>(I_amp / corr); // = G*Itrue
|
|
r.sigma = static_cast<float>(sigma_amp / corr);
|
|
r.image_scale_corr = static_cast<float>(1.0 / data.scale_factor); // G only
|
|
} else {
|
|
r.I = static_cast<float>(I_amp);
|
|
r.sigma = static_cast<float>(sigma_amp);
|
|
r.image_scale_corr = NAN;
|
|
}
|
|
} else {
|
|
r.I = 0.0f;
|
|
r.sigma = NAN;
|
|
r.bkg = 0.0f;
|
|
r.observed = false;
|
|
}
|
|
data.reflections.push_back(r);
|
|
}
|
|
|
|
// ---- Per-image CC vs reference (the half/ref correlation diagnostic) -------
|
|
// Pearson CC of the scaled estimate (r.I * image_scale_corr = Itrue_est)
|
|
// against the reference intensities, over the matched reflections.
|
|
{
|
|
double sx = 0, sy = 0, sxx = 0, syy = 0, sxy = 0;
|
|
size_t cn = 0;
|
|
for (const auto &r : data.reflections) {
|
|
if (!r.observed || !std::isfinite(r.I) || !std::isfinite(r.image_scale_corr))
|
|
continue;
|
|
const auto it = reference_data.find(hkl_key_generator(r));
|
|
if (it == reference_data.end())
|
|
continue;
|
|
const double x = static_cast<double>(r.I) * static_cast<double>(r.image_scale_corr);
|
|
const double y = it->second;
|
|
if (!std::isfinite(x) || !std::isfinite(y))
|
|
continue;
|
|
sx += x; sy += y; sxx += x * x; syy += y * y; sxy += x * y; ++cn;
|
|
}
|
|
data.cc = NAN;
|
|
data.cc_n = static_cast<int64_t>(cn);
|
|
if (cn >= 3) {
|
|
const double nd = static_cast<double>(cn);
|
|
const double cov = sxy - sx * sy / nd;
|
|
const double vx = sxx - sx * sx / nd;
|
|
const double vy = syy - sy * sy / nd;
|
|
if (vx > 0.0 && vy > 0.0)
|
|
data.cc = cov / std::sqrt(vx * vy);
|
|
}
|
|
}
|
|
}
|
|
|
|
std::vector<float> PixelRefine::PredictImage(const AzimuthalIntegrationProfile &profile,
|
|
BraggPrediction &prediction,
|
|
const PixelRefineData &data,
|
|
bool include_background) const {
|
|
std::vector<float> img(xpixel * ypixel, 0.0f);
|
|
|
|
const double lambda = data.geom.GetWavelength_A();
|
|
const double pixel_size = data.geom.GetPixelSize_mm();
|
|
const auto azim_result = profile.GetResult();
|
|
const auto &pixel_to_bin = mapping.GetPixelToBin();
|
|
const auto &corrections = mapping.Corrections();
|
|
const int total_bin_count = static_cast<int>(azim_result.size());
|
|
const double angle_rad = data.angle_deg * M_PI / 180.0;
|
|
const int radius = data.shoebox_radius;
|
|
const double bw = data.bandwidth;
|
|
|
|
auto recip_area = [&](double x, double y) -> double {
|
|
const Coord qx = data.geom.DetectorToRecip(x + 0.5, y) - data.geom.DetectorToRecip(x - 0.5, y);
|
|
const Coord qy = data.geom.DetectorToRecip(x, y + 0.5) - data.geom.DetectorToRecip(x, y - 0.5);
|
|
return (qx % qy).Length();
|
|
};
|
|
auto bandwidth_radial_sq = [&](double d) -> double {
|
|
if (bw <= 0.0 || d <= 0.0)
|
|
return 0.0;
|
|
const double bl = bw * lambda;
|
|
return bl * bl / (2.0 * d * d * d * d);
|
|
};
|
|
|
|
// The model works in solid-angle/polarization-corrected units (as in Run,
|
|
// where Iobs = raw * correction). Map back to raw detector units (/ correction)
|
|
// so the predicted image overlays directly on the original image.
|
|
auto to_raw = [&](size_t npixel, double corrected) -> float {
|
|
const double corr = corrections[npixel];
|
|
return (corr > 0.0) ? static_cast<float>(corrected / corr) : 0.0f;
|
|
};
|
|
|
|
// Background base layer (per-pixel azimuthal mean), full-frame pass.
|
|
if (include_background) {
|
|
for (size_t p = 0; p < img.size(); ++p) {
|
|
const int bin = pixel_to_bin[p];
|
|
if (bin >= 0 && bin < total_bin_count)
|
|
img[p] = to_raw(p, azim_result[bin]);
|
|
}
|
|
}
|
|
|
|
double beam[2], dist_mm, detector_rot[2], rot_vec[3];
|
|
double latt_vec0[3], latt_vec1[3], latt_vec2[3];
|
|
BuildParameterBlocks(data, beam, dist_mm, detector_rot, rot_vec, latt_vec0, latt_vec1, latt_vec2);
|
|
|
|
DiffractionExperiment exp_iter = experiment;
|
|
exp_iter.BeamX_pxl(data.geom.GetBeamX_pxl())
|
|
.BeamY_pxl(data.geom.GetBeamY_pxl())
|
|
.DetectorDistance_mm(data.geom.GetDetectorDistance_mm())
|
|
.PoniRot1_rad(data.geom.GetPoniRot1_rad())
|
|
.PoniRot2_rad(data.geom.GetPoniRot2_rad());
|
|
|
|
const BraggPredictionSettings settings_prediction{
|
|
.high_res_A = experiment.GetBraggIntegrationSettings().GetDMinLimit_A(),
|
|
.max_hkl = 100,
|
|
.centering = data.centering,
|
|
.bandwidth_sigma = static_cast<float>(data.bandwidth) // relative Δλ/λ sigma
|
|
};
|
|
const int nrefl = prediction.Calc(exp_iter, data.latt, settings_prediction);
|
|
const auto &predicted = prediction.GetReflections();
|
|
|
|
for (int ri = 0; ri < nrefl; ++ri) {
|
|
const auto &refl = predicted[ri];
|
|
const auto it = reference_data.find(hkl_key_generator(refl));
|
|
if (it == reference_data.end())
|
|
continue;
|
|
|
|
const double Itrue = it->second;
|
|
const double R_bw_sq = bandwidth_radial_sq(refl.d);
|
|
|
|
const int min_y = std::max<int>(refl.predicted_y - radius, 0);
|
|
const int max_y = std::min<int>(refl.predicted_y + radius, ypixel - 1);
|
|
const int min_x = std::max<int>(refl.predicted_x - radius, 0);
|
|
const int max_x = std::min<int>(refl.predicted_x + radius, xpixel - 1);
|
|
|
|
for (int y = min_y; y <= max_y; ++y) {
|
|
for (int x = min_x; x <= max_x; ++x) {
|
|
const size_t npixel = xpixel * y + x;
|
|
|
|
// Pure Bragg signal: Ibkg = 0 so Model() returns signal only; the
|
|
// background is already laid down above. Same code path as Run.
|
|
PixelObs obs{
|
|
.x = static_cast<double>(x),
|
|
.y = static_cast<double>(y),
|
|
.Iobs = 0.0,
|
|
.Ibkg = 0.0,
|
|
.weight = 1.0,
|
|
.A_recip = recip_area(x, y),
|
|
.angle_rad = angle_rad
|
|
};
|
|
PixelResidual pr(obs, Itrue, lambda, pixel_size,
|
|
refl.h, refl.k, refl.l, R_bw_sq, data.crystal_system);
|
|
|
|
double signal = 0.0;
|
|
if (pr.Model(beam, &dist_mm, detector_rot, rot_vec,
|
|
latt_vec0, latt_vec1, latt_vec2,
|
|
&data.scale_factor, &data.B_factor, data.R, signal))
|
|
img[npixel] += to_raw(npixel, signal);
|
|
}
|
|
}
|
|
}
|
|
|
|
return img;
|
|
}
|
|
|
|
template<class T>
|
|
std::vector<float> PixelRefine::ChiSquaredImage(const T *image,
|
|
const AzimuthalIntegrationProfile &profile,
|
|
BraggPrediction &prediction,
|
|
const PixelRefineData &data) const {
|
|
std::vector<float> img(xpixel * ypixel, 0.0f);
|
|
|
|
const double lambda = data.geom.GetWavelength_A();
|
|
const double pixel_size = data.geom.GetPixelSize_mm();
|
|
const auto azim_result = profile.GetResult();
|
|
const auto azim_std = profile.GetStd();
|
|
const auto &pixel_to_bin = mapping.GetPixelToBin();
|
|
const auto &corrections = mapping.Corrections();
|
|
const int total_bin_count = static_cast<int>(azim_result.size());
|
|
const double angle_rad = data.angle_deg * M_PI / 180.0;
|
|
const int radius = data.shoebox_radius;
|
|
const double bw = data.bandwidth;
|
|
|
|
auto recip_area = [&](double x, double y) -> double {
|
|
const Coord qx = data.geom.DetectorToRecip(x + 0.5, y) - data.geom.DetectorToRecip(x - 0.5, y);
|
|
const Coord qy = data.geom.DetectorToRecip(x, y + 0.5) - data.geom.DetectorToRecip(x, y - 0.5);
|
|
return (qx % qy).Length();
|
|
};
|
|
auto bandwidth_radial_sq = [&](double d) -> double {
|
|
if (bw <= 0.0 || d <= 0.0)
|
|
return 0.0;
|
|
const double bl = bw * lambda;
|
|
return bl * bl / (2.0 * d * d * d * d);
|
|
};
|
|
|
|
double beam[2], dist_mm, detector_rot[2], rot_vec[3];
|
|
double latt_vec0[3], latt_vec1[3], latt_vec2[3];
|
|
BuildParameterBlocks(data, beam, dist_mm, detector_rot, rot_vec, latt_vec0, latt_vec1, latt_vec2);
|
|
|
|
DiffractionExperiment exp_iter = experiment;
|
|
exp_iter.BeamX_pxl(data.geom.GetBeamX_pxl())
|
|
.BeamY_pxl(data.geom.GetBeamY_pxl())
|
|
.DetectorDistance_mm(data.geom.GetDetectorDistance_mm())
|
|
.PoniRot1_rad(data.geom.GetPoniRot1_rad())
|
|
.PoniRot2_rad(data.geom.GetPoniRot2_rad());
|
|
|
|
const BraggPredictionSettings settings_prediction{
|
|
.high_res_A = experiment.GetBraggIntegrationSettings().GetDMinLimit_A(),
|
|
.max_hkl = 100,
|
|
.centering = data.centering,
|
|
.bandwidth_sigma = static_cast<float>(data.bandwidth)
|
|
};
|
|
const int nrefl = prediction.Calc(exp_iter, data.latt, settings_prediction);
|
|
const auto &predicted = prediction.GetReflections();
|
|
|
|
for (int ri = 0; ri < nrefl; ++ri) {
|
|
const auto &refl = predicted[ri];
|
|
const auto it = reference_data.find(hkl_key_generator(refl));
|
|
if (it == reference_data.end())
|
|
continue;
|
|
|
|
const double Itrue = it->second;
|
|
const double R_bw_sq = bandwidth_radial_sq(refl.d);
|
|
|
|
const int min_y = std::max<int>(refl.predicted_y - radius, 0);
|
|
const int max_y = std::min<int>(refl.predicted_y + radius, ypixel - 1);
|
|
const int min_x = std::max<int>(refl.predicted_x - radius, 0);
|
|
const int max_x = std::min<int>(refl.predicted_x + radius, xpixel - 1);
|
|
|
|
for (int y = min_y; y <= max_y; ++y) {
|
|
for (int x = min_x; x <= max_x; ++x) {
|
|
const size_t npixel = xpixel * y + x;
|
|
const int azim_bin = pixel_to_bin[npixel];
|
|
|
|
// Same gating as Run(): only pixels that actually enter the fit.
|
|
if (azim_bin >= total_bin_count)
|
|
continue;
|
|
if (image[npixel] == std::numeric_limits<T>::max())
|
|
continue;
|
|
if (std::is_signed_v<T> && (image[npixel] == std::numeric_limits<T>::min()))
|
|
continue;
|
|
|
|
const double correction = corrections[npixel];
|
|
const double Ibkg = azim_result[azim_bin];
|
|
const double Ibkg_sigma = azim_std[azim_bin];
|
|
const double raw = static_cast<double>(image[npixel]);
|
|
const double Iobs = raw * correction;
|
|
|
|
double var = correction * std::max(Iobs, 0.0) + Ibkg_sigma * Ibkg_sigma;
|
|
if (!(var > 1.0))
|
|
var = 1.0;
|
|
const double weight = 1.0 / std::sqrt(var);
|
|
|
|
PixelObs obs{
|
|
.x = static_cast<double>(x),
|
|
.y = static_cast<double>(y),
|
|
.Iobs = Iobs,
|
|
.Ibkg = Ibkg,
|
|
.weight = weight,
|
|
.A_recip = recip_area(x, y),
|
|
.angle_rad = angle_rad
|
|
};
|
|
PixelResidual pr(obs, Itrue, lambda, pixel_size,
|
|
refl.h, refl.k, refl.l, R_bw_sq, data.crystal_system);
|
|
|
|
double Ipred = 0.0;
|
|
if (pr.Model(beam, &dist_mm, detector_rot, rot_vec,
|
|
latt_vec0, latt_vec1, latt_vec2,
|
|
&data.scale_factor, &data.B_factor, data.R, Ipred)) {
|
|
// residual_i = (I_pred - I_obs) * weight (== Ceres residual);
|
|
// its square is this pixel's contribution to the cost.
|
|
const double rw = (Ipred - Iobs) * weight;
|
|
img[npixel] += static_cast<float>(rw * rw);
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
return img;
|
|
}
|
|
|
|
// Explicit instantiations for the supported (uncompressed) image pixel types.
|
|
template void PixelRefine::Run<int8_t>(const int8_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, PixelRefineData &);
|
|
template void PixelRefine::Run<int16_t>(const int16_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, PixelRefineData &);
|
|
template void PixelRefine::Run<int32_t>(const int32_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, PixelRefineData &);
|
|
template void PixelRefine::Run<uint8_t>(const uint8_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, PixelRefineData &);
|
|
template void PixelRefine::Run<uint16_t>(const uint16_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, PixelRefineData &);
|
|
template void PixelRefine::Run<uint32_t>(const uint32_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, PixelRefineData &);
|
|
|
|
template std::vector<float> PixelRefine::ChiSquaredImage<int8_t>(const int8_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, const PixelRefineData &) const;
|
|
template std::vector<float> PixelRefine::ChiSquaredImage<int16_t>(const int16_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, const PixelRefineData &) const;
|
|
template std::vector<float> PixelRefine::ChiSquaredImage<int32_t>(const int32_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, const PixelRefineData &) const;
|
|
template std::vector<float> PixelRefine::ChiSquaredImage<uint8_t>(const uint8_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, const PixelRefineData &) const;
|
|
template std::vector<float> PixelRefine::ChiSquaredImage<uint16_t>(const uint16_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, const PixelRefineData &) const;
|
|
template std::vector<float> PixelRefine::ChiSquaredImage<uint32_t>(const uint32_t *, const AzimuthalIntegrationProfile &, BraggPrediction &, const PixelRefineData &) const;
|