Files
Jungfraujoch/image_analysis/RotationIndexer.cpp
Filip Leonarski 07fe4dd3bb
All checks were successful
Build Packages / build:rpm (rocky9_nocuda) (push) Successful in 11m23s
Build Packages / build:rpm (ubuntu2204_nocuda) (push) Successful in 10m32s
Build Packages / build:rpm (ubuntu2404_nocuda) (push) Successful in 9m15s
Build Packages / Generate python client (push) Successful in 19s
Build Packages / Build documentation (push) Successful in 49s
Build Packages / Create release (push) Has been skipped
Build Packages / build:rpm (rocky8_sls9) (push) Successful in 9m13s
Build Packages / build:rpm (rocky8) (push) Successful in 9m10s
Build Packages / build:rpm (rocky9) (push) Successful in 9m58s
Build Packages / build:rpm (ubuntu2204) (push) Successful in 8m52s
Build Packages / build:rpm (ubuntu2404) (push) Successful in 8m42s
Build Packages / Unit tests (push) Successful in 1h12m44s
Build Packages / build:rpm (rocky8_nocuda) (push) Successful in 11m30s
v1.0.0-rc.124 (#31)
This is an UNSTABLE release. This version significantly rewrites code to predict reflection position and integrate them,
especially in case of rotation crystallography. If things go wrong with analysis, it is better to revert to 1.0.0-rc.123.

* jfjoch_broker: Improve refection position prediction and Bragg integration code.
* jfjoch_broker: Align with XDS way of calculating Lorentz correction and general notation.
* jfjoch_writer: Fix saving mosaicity properly in HDF5 file.
* jfjoch_viewer: Introduce high-dynamic range mode for images
* jfjoch_viewer: Ctrl+mouse wheel has exponential change in foreground (+/-15%)
* jfjoch_viewer: Zoom-in numbers have better readability

Reviewed-on: #31
Co-authored-by: Filip Leonarski <filip.leonarski@psi.ch>
Co-committed-by: Filip Leonarski <filip.leonarski@psi.ch>
2026-02-01 13:29:33 +01:00

153 lines
5.0 KiB
C++

// SPDX-FileCopyrightText: 2025 Filip Leonarski, Paul Scherrer Institute <filip.leonarski@psi.ch>
// SPDX-License-Identifier: GPL-3.0-only
#include "RotationIndexer.h"
#include "geom_refinement/XtalOptimizer.h"
#include "indexing/FFTIndexer.h"
#include "lattice_search/LatticeSearch.h"
#include <iostream>
RotationIndexer::RotationIndexer(const DiffractionExperiment &x, IndexerThreadPool &indexer)
: experiment(x),
index_ice_rings(x.GetIndexingSettings().GetIndexIceRings()),
axis_(x.GetGoniometer()),
geom_(x.GetDiffractionGeometry()),
updated_geom_(geom_),
indexer_(indexer) {
if (axis_) {
angle_norm_deg = std::fabs(axis_->GetIncrement_deg());
if (angle_norm_deg < 1e-6) {
// Guard against rotation close to zero
axis_ = std::nullopt;
} else {
if (x.GetImageNum() < min_images_for_indexing) {
// For short measurements - only indexing at the end
first_image_to_try_indexing = INT64_MAX;
image_stride = 1;
} else {
first_image_to_try_indexing = std::max<int64_t>(min_images_for_indexing,
x.GetIndexingSettings().GetRotationIndexingMinAngularRange_deg() / angle_norm_deg);
image_stride = std::ceil(x.GetIndexingSettings().GetRotationIndexingAngularStride_deg() / angle_norm_deg);
if (image_stride == 0)
image_stride = 1;
}
}
}
}
void RotationIndexer::TryIndex() {
indexing_tried = true;
// Index
std::vector<SpotToSave> v_sel;
std::vector<Coord> coords_sel;
if (coords_.size() > max_spots) {
// Indices into v_ / coords_
std::vector<std::size_t> idx(coords_.size());
std::iota(idx.begin(), idx.end(), std::size_t{0});
// Sort indices by descending intensity
std::partial_sort(
idx.begin(),
idx.begin() + max_spots,
idx.end(),
[this](std::size_t a, std::size_t b) {
return v_[a].intensity > v_[b].intensity;
}
);
v_sel.reserve(max_spots);
coords_sel.reserve(max_spots);
for (std::size_t i = 0; i < max_spots; ++i) {
const std::size_t k = idx[i];
v_sel.emplace_back(v_[k]);
coords_sel.emplace_back(coords_[k]);
}
} else {
v_sel = v_;
coords_sel = coords_;
}
auto indexer_result = indexer_.Run(experiment, coords_sel).get();
if (!indexer_result.lattice.empty()) {
// Find lattice type
search_result_ = LatticeSearch(indexer_result.lattice[0]);
// Run refinement
DiffractionExperiment experiment_copy(experiment);
XtalOptimizerData data{
.geom = experiment_copy.GetDiffractionGeometry(),
.latt = search_result_.conventional,
.crystal_system = search_result_.system,
.min_spots = experiment.GetIndexingSettings().GetViableCellMinSpots(),
.refine_beam_center = true,
.refine_distance_mm = false,
.refine_detector_angles = true,
.refine_rotation_axis = true,
.index_ice_rings = experiment.GetIndexingSettings().GetIndexIceRings(),
.axis = axis_
};
if (data.crystal_system == gemmi::CrystalSystem::Trigonal)
data.crystal_system = gemmi::CrystalSystem::Hexagonal;
if (XtalOptimizer(data, v_sel)) {
indexed_lattice = data.latt;
updated_geom_ = data.geom;
axis_ = data.axis;
}
}
}
std::optional<RotationIndexerResult> RotationIndexer::ProcessImage(int64_t image, const std::vector<SpotToSave> &spots) {
std::unique_lock ul(m);
// For non-rotation just ignore the whole procedure
if (!axis_)
return {};
const float angle_deg = axis_->GetAngle_deg(image) + axis_->GetWedge_deg() / 2.0f;
const auto rot = axis_->GetTransformationAngle(angle_deg);
if (!indexing_tried && image >= last_accumulated_image + image_stride) {
v_.reserve(v_.size() + spots.size());
coords_.reserve(coords_.size() + spots.size());
for (const auto &s: spots) {
if (index_ice_rings || !s.ice_ring) {
v_.emplace_back(s);
coords_.emplace_back(rot * s.ReciprocalCoord(geom_));
}
}
accumulated_images++;
last_accumulated_image = image;
if (accumulated_images >= min_images_for_indexing && image >= first_image_to_try_indexing)
TryIndex();
}
if (!indexed_lattice)
return {};
return RotationIndexerResult{
.lattice = indexed_lattice.value(),
.search_result = search_result_,
.geom = updated_geom_,
.axis = axis_
};
}
std::optional<RotationIndexerResult> RotationIndexer::GetLattice() {
if (!indexed_lattice)
return {};
return RotationIndexerResult{
.lattice = indexed_lattice.value(),
.search_result = search_result_,
.geom = updated_geom_,
.axis = axis_
};
}