Code Reference¶
The following is not an exhaustive code reference, but a list of the C++ objects exposed via the Python API. These objects can be used to perform a complete data reduction via a Python script or (for example) a Jupyter notebook.
- group python_api
Typedefs
-
using RankedSolution = std::pair<sptrUnitCell, double>¶
-
using SolutionList = std::vector<RankedSolution>¶
-
using EquivalencePair = std::pair<int, int>¶
-
using EquivalenceList = std::vector<EquivalencePair>¶
-
typedef std::shared_ptr<ProgressHandler> sptrProgressHandler¶
-
using InstrumentStateList = std::vector<InstrumentState>¶
-
typedef std::shared_ptr<Peak3D> sptrPeak3D¶
-
typedef std::vector<sptrPeak3D> PeakList¶
-
using sptrShapeModel = std::shared_ptr<ShapeModel>¶
Enums
-
enum class RejectionFlag¶
Specifies reason why a peak was rejected.
Values:
-
enumerator NotRejected¶
-
enumerator Masked¶
-
enumerator OutsideThreshold¶
-
enumerator OutsideFrames¶
-
enumerator OutsideDetector¶
-
enumerator IntegrationFailure¶
-
enumerator TooFewPoints¶
-
enumerator NoNeighbours¶
-
enumerator TooFewNeighbours¶
-
enumerator NoUnitCell¶
-
enumerator NoDataSet¶
-
enumerator InvalidRegion¶
-
enumerator InterpolationFailure¶
-
enumerator InvalidShape¶
-
enumerator InvalidSigma¶
-
enumerator InvalidBkgSigma¶
-
enumerator SaturatedPixel¶
-
enumerator OverlappingBkg¶
-
enumerator OverlappingPeak¶
-
enumerator InvalidCentroid¶
-
enumerator InvalidCovariance¶
-
enumerator CentreOutOfBounds¶
-
enumerator BadIntegrationFit¶
-
enumerator NoShapeModel¶
-
enumerator NoISigmaMinimum¶
-
enumerator TooWide¶
-
enumerator BadGaussianFit¶
-
enumerator PredictionUpdateFailure¶
-
enumerator ManuallyRejected¶
-
enumerator OutsideIndexingTol¶
-
enumerator Outlier¶
-
enumerator Extinct¶
-
enumerator Count¶
-
enumerator NotRejected¶
-
enum class PeakCollectionType¶
Type of peak collection.
Values:
-
enumerator FOUND¶
-
enumerator PREDICTED¶
-
enumerator INDEXING¶
-
enumerator FOUND¶
Functions
-
void log(const Level &level) const¶
-
Refiner(UnitCellHandler *cell_handler)¶
-
void reconstructBatches(std::vector<Peak3D*> peaks)¶
Rebuild old batches if refinement failed.
-
void refineUB()¶
Sets the lattice B matrix to be refined.
-
void refineDetectorOffset()¶
Sets detector offsets in the given list of instrument states to be refined.
-
void refineSamplePosition()¶
Sets the sample position in the given list of instrument states to be refined.
-
void refineSampleOrientation()¶
Sets the sample orientation in the given list of instrument states to be refined.
-
void refineKi()¶
Sets the source ki in the given list of instrument states to be refined.
-
bool refine(sptrDataSet data, const std::vector<Peak3D*> peaks, sptrUnitCell cell = nullptr)¶
Perform the refinement with the maximum number of iterations as given. N.B. the four previous functions set the number of free parameters and at least one must be run before refine
-
int updatePredictions(std::vector<Peak3D*> peaks)¶
Updates the centers of predicted peaks, after refinement. Returns the number of remaining peaks (some peaks may be rejected with flag PredictionUpdateFailure).
-
const std::vector<RefinementBatch> &batches() const¶
Returns the individual peak/frame batches used during refinement.
-
const UnitCell *unrefinedCell() const¶
Return the unrefined cell.
-
const InstrumentStateList *refinedStates() const¶
Return the refined states.
-
const InstrumentStateList *unrefinedStates() const¶
Return the unrefined states.
-
int nframes() const¶
Return number of frames.
-
bool firstRefine() const¶
Check if this is the first refinement.
-
void logChange()¶
Write the initial and final cells to the log.
-
RefinerParameters *parameters()¶
get a pointer to the parameters
-
void setParameters(const RefinerParameters ¶ms)¶
set the parameters
-
void assignPredictedCells(std::vector<Peak3D*> predicted_peaks)¶
Assign batch cells to predicted peaks.
-
void makeBatches(InstrumentStateList &states, const std::vector<Peak3D*> &peaks, sptrUnitCell cell)¶
Generate batches of peaks per frame range with the given peak list.
-
sptrUnitCell _getUnitCell(const std::vector<Peak3D*> peaks_subset)¶
Determine which unit cell to use in a batch.
-
ExperimentYAML(const std::string &filename)¶
-
ExperimentYAML(const YAML::Node &node)¶
-
void grabDataReaderParameters(DataReaderParameters *params) const¶
Get parameters from data reader.
-
void setDataReaderParameters(DataReaderParameters *params)¶
Set data reader parameters.
-
void grabIntegrationParameters(IntegrationParameters *params)¶
Get parameters from integrator.
-
void setIntegrationParameters(IntegrationParameters *params)¶
Set integrator parameters.
-
void grabPeakFinderParameters(PeakFinderParameters *params)¶
Get parameters from integrator.
-
void setPeakFinderParameters(PeakFinderParameters *params)¶
Set peak finder parameters.
-
void grabAutoindexerParameters(IndexerParameters *params)¶
Get parameters from indexer.
-
void setAutoindexerParameters(IndexerParameters *params)¶
Set indexer parameters.
-
void grabShapeParameters(ShapeModelParameters *params)¶
Get parameters from shape model.
-
void setShapeParameters(ShapeModelParameters *params)¶
Set shape model parameters.
-
void grabPredictorParameters(PredictionParameters *params)¶
Get parameters from predictor.
-
void setPredictorParameters(PredictionParameters *params)¶
Set shape model parameters.
-
void grabMergeParameters(MergeParameters *params)¶
Get parameters from merger.
-
void setMergeParameters(MergeParameters *params)¶
Set merger parameters.
-
void writeFile(const std::string &filename)¶
Write the yml file.
-
template<typename T>
inline T getNode(const YAML::Node &node, const std::string &key) const¶ Read an arbitrary type node, check whether it exists.
-
std::ostream &operator<<(std::ostream &os, const Blob3D &b)¶
Variables
-
bool refine_ub = true¶
Refine unit cell.
-
bool refine_ki = true¶
Refine incident wavevector.
-
bool refine_sample_position = true¶
Refine sample position.
-
bool refine_sample_orientation = true¶
Refine sample orientation.
-
bool refine_detector_offset = true¶
Refine detector offset.
-
int nbatches = 10¶
Number of refinement batches.
-
ResidualType residual_type = ResidualType::QSpace¶
Type of residual used in minimisation.
-
unsigned int max_iter = 1000¶
Maximum number of iteration for minimisation.
-
bool use_batch_cells = false¶
Use unit cells assigned to batches in further refinement.
-
bool set_unit_cell = true¶
Assign unit cell to peaks when making batches.
-
UnitCellHandler *_cell_handler¶
-
UnitCell _unrefined_cell¶
-
InstrumentStateList _unrefined_states¶
-
sptrUnitCell _cell¶
-
std::vector<UnitCell*> _batch_cells¶
-
std::vector<RefinementBatch> _batches¶
-
int _nframes¶
-
bool _first_refine = true¶
-
InstrumentStateList *_states¶
-
std::unique_ptr<RefinerParameters> _params¶
-
static constexpr double _eps_norm = 50.0¶
-
std::vector<sptrUnitCell> _tmp_vec¶
-
std::vector<Peak3D*> _peaks¶
-
YAML::Node _node¶
The root yml node.
-
struct IndexerParameters¶
- #include <AutoIndexer.h>
Parameters for indexing a peak collection.
Public Members
-
double maxdim = 200.0¶
Expected maximum dimension of the unit cell, in angstroms.
-
int nSolutions = 10¶
Number of candidate lattice vectors to use to search for a unit cell.
-
int nVertices = 10000¶
Number of points on the unit sphere to sample when looking for lattice directions
-
int subdiv = 30¶
Parameter used to control the number of histogram bins to use during FFT.
-
double indexingTolerance = 0.2¶
Tolerance used to determine if a peak is indexed by the unit cell.
-
double niggliTolerance = 1e-3¶
Tolerance to use when reducing to a Niggli cell.
-
double gruberTolerance = 4e-2¶
Tolerance to use during Gruber reduction.
-
bool niggliReduction = false¶
Use only Niggli reduction (i.e., no subsequent Gruber reduction) if set true
-
double minUnitCellVolume = 20.0¶
Lower bound of expected unit cell volume, used to reject candidate unit cells
-
double unitCellEquivalenceTolerance = 0.05¶
Tolerance value used to decide of two unit cells are equivalent.
-
double solutionCutoff = 10.0¶
Indexing quality cutoff (value indicates percentage of peaks indexed by cell)
-
double frequencyTolerance = 0.7¶
Frequency threshold: FFT peaks with 0.7 the value of the zero frequency peak are discarded
-
double first_frame = 0.0¶
Frame from which to start autoindexing.
-
double last_frame = 10.0¶
Last frame of autoindexing set.
-
double d_min = 1.5¶
Minimum detector range.
-
double d_max = 50.0¶
Maximum detector range.
-
double strength_min = 1.0¶
Minimum peak strength.
-
double strength_max = 1.0e7¶
Maximum peak strength.
-
bool peaks_integrated = false¶
Whether the peaks are integrated.
-
double maxdim = 200.0¶
-
class AutoIndexer¶
- #include <AutoIndexer.h>
Class to automatically index a set of peaks.
This class uses 1D Fourier transform algorithm to generate a list of candidate unit cells (I. Steller et al., J. Appl. Cryst., 30:1036, 1997. doi:10.1107/S0021889897008777).
A set of random q vectors is generated on the unit sphere, and project the peak q vectors onto each direction; if the direction corresponds to a lattice plane, the Fourier transform will be strongly periodic, and the first non-zero frequency corresponds to the lattice plane. Each peak is assigned a set of Miller indices and the unit cell determined from these.
Public Functions
-
AutoIndexer()¶
Constructor.
-
IndexerParameters *parameters()¶
Return the autoindexing parameters.
-
bool autoIndex(const std::vector<Peak3D*> &peaks, sptrDataSet data, const InstrumentState *state = nullptr, bool filter = true)¶
Perform the autoindexing, possibly for a single frame only.
-
bool autoIndex(PeakCollection *peaks, sptrDataSet data, const InstrumentState *state = nullptr, bool filter = true)¶
Autoindex by passing a peak collection (avoid SWIG memory leak)
-
const std::vector<RankedSolution> &solutions() const¶
Return a list of the best solutions ordered by percentage of successfully indexed peaks.
Set the progress handler.
-
inline void unsetHandler()¶
Remove the progress handler.
-
std::string solutionsToString() const¶
Dump SolutionList as a string.
-
void acceptSolution(const sptrUnitCell solution, const std::vector<ohkl::Peak3D*> &peaks)¶
Set solution to be unit cell for given peak list.
-
void acceptSolution(const sptrUnitCell solution, PeakCollection *peaks)¶
Set solution to be unit cell for given peak collection.
-
sptrUnitCell goodSolution(const UnitCell *reference_cell, double length_tol, double angle_tol)¶
Check if list of solutions contains reference unit cell. If it does, return a pointer, otherwise return nullptr
-
std::vector<RankedSolution> filterSolutions() const¶
Filter out solutions that do not meet criteria.
-
AutoIndexer()¶
-
struct RefinerParameters¶
- #include <Refiner.h>
Parameters for refinement.
Public Members
-
bool refine_ub = true¶
Refine unit cell.
-
bool refine_ki = true¶
Refine incident wavevector.
-
bool refine_sample_position = true¶
Refine sample position.
-
bool refine_sample_orientation = true¶
Refine sample orientation.
-
bool refine_detector_offset = true¶
Refine detector offset.
-
int nbatches = 10¶
Number of refinement batches.
-
ResidualType residual_type = ResidualType::QSpace¶
Type of residual used in minimisation.
-
unsigned int max_iter = 1000¶
Maximum number of iteration for minimisation.
-
bool use_batch_cells = false¶
Use unit cells assigned to batches in further refinement.
-
bool set_unit_cell = true¶
Assign unit cell to peaks when making batches.
-
bool refine_ub = true¶
-
class Refiner¶
- #include <Refiner.h>
Refine unit cell and instrument states.
Use GSL non-linear least squares minimisation to refine instrument states (sample position and orientation, incident wavevector, and detector position) and unit cell parameters. The residual is computed as the difference between a q-vector computed from the unit cell (used to predict peaks) and a q-vector from transforming an observed reflection into reciprocal space. An alternative residual (useful for refining unit cells and incident wavevectors) is computed using the real space difference between a predicted spot and its associated observed spot.
Public Functions
-
void refineUB()¶
Sets the lattice B matrix to be refined.
-
void refineDetectorOffset()¶
Sets detector offsets in the given list of instrument states to be refined.
-
void refineSamplePosition()¶
Sets the sample position in the given list of instrument states to be refined.
-
void refineSampleOrientation()¶
Sets the sample orientation in the given list of instrument states to be refined.
-
void refineKi()¶
Sets the source ki in the given list of instrument states to be refined.
-
bool refine(sptrDataSet data, const std::vector<Peak3D*> peaks, sptrUnitCell cell = nullptr)¶
Perform the refinement with the maximum number of iterations as given. N.B. the four previous functions set the number of free parameters and at least one must be run before refine
-
int updatePredictions(std::vector<Peak3D*> peaks)¶
Updates the centers of predicted peaks, after refinement. Returns the number of remaining peaks (some peaks may be rejected with flag PredictionUpdateFailure).
-
const std::vector<RefinementBatch> &batches() const¶
Returns the individual peak/frame batches used during refinement.
-
const UnitCell *unrefinedCell() const¶
Return the unrefined cell.
-
const InstrumentStateList *refinedStates() const¶
Return the refined states.
-
const InstrumentStateList *unrefinedStates() const¶
Return the unrefined states.
-
int nframes() const¶
Return number of frames.
-
bool firstRefine() const¶
Check if this is the first refinement.
-
void logChange()¶
Write the initial and final cells to the log.
-
RefinerParameters *parameters()¶
get a pointer to the parameters
-
void setParameters(const RefinerParameters ¶ms)¶
set the parameters
-
void refineUB()¶
-
class DataSet¶
- #include <DataSet.h>
Class used to manage loading of detector images and metadata from disk.
Note that this class does not contain any actual data, just a few metadata and a link to the data reader.
Subclassed by ohkl::SingleFrame
Public Functions
-
std::size_t nFrames() const¶
Number of detector image frames acquired so far.
-
std::size_t nRows() const¶
The number of rows in each detector image.
-
std::size_t nCols() const¶
The number of columns in each detector image.
-
void addBoxMask(const AABB &aabb)¶
Add a rectangular mask (swig-compatible)
-
void addEllipseMask(const AABB &aabb)¶
Add an elliptical mask (swig-compatible)
-
virtual Eigen::MatrixXi frame(const std::size_t idx) const¶
Read a single frame.
-
virtual Eigen::MatrixXd transformedFrame(std::size_t idx) const¶
Returns frame after transforming to account for detector gain and baseline.
-
virtual Eigen::MatrixXd gradientFrame(std::size_t idx, GradientKernel kernel, bool realspace = true) const¶
Return per-pixel magnitude of gradient of a given frame.
-
void open()¶
Gets the file handle.
-
void close()¶
Close file and release handle.
-
ReciprocalVector computeQ(const DetectorEvent &ev) const¶
Returns the sample-space q vector corresponding to a detector event.
-
const IDataReader *reader() const¶
Returns the data reader used to set this dataset.
-
const Diffractometer *diffractometer() const¶
Returns the diffractometer associated to this dataset.
-
Diffractometer *diffractometer()¶
Returns the diffractometer associated to this dataset.
-
Detector &detector()¶
Returns the detector associated to this dataset.
-
const Detector &detector() const¶
Returns the detector associated to this dataset.
-
const ohkl::MetaData &metadata() const¶
Returns a const reference to the MetaData container.
-
ohkl::MetaData &metadata()¶
Returns a reference to the MetaData container.
-
void addDataFile(const std::string &filename, const DataFormat &fmt)¶
Add a data file for reading data. Reading frames will be done only upon request.
-
void setImageReaderParameters(const DataReaderParameters ¶ms)¶
Set the parameters for the image (raw or tiff) reader.
-
void addFrame(const std::string &filename, const DataFormat &format)¶
Add an image from a sequence.
-
void finishRead()¶
Finish reading procedure (must be called before using the data stored in the DataSet).
-
double wavelength() const¶
Query the wavelength stored in the metadata.
-
BitDepth bitDepth() const¶
Query image bit depth stored in metadata.
-
void setInstrumentStates(InstrumentStateSet *states)¶
Get the initial instrument states.
-
InstrumentStateList &instrumentStates()¶
Return instrument state list.
-
void adjustDirectBeam(double x_offset, double y_offset)¶
Adjust the direct beam (mainly for scripting purposes)
-
void initHistograms(std::size_t nbins)¶
Initialise intensity histograms.
-
void getFrameIntensityHistogram(std::size_t index)¶
Generate frame intensity histogram.
-
void clearHistograms()¶
Free histogram memeory.
-
inline size_t getNumberHistograms()¶
Get number of available histograms.
-
gsl_histogram *getHistogram(int index)¶
accessing created histograms
-
gsl_histogram *getTotalHistogram()¶
accessing Total histogram
-
bool hasMasks()¶
returns a booleans whether masks have been created or not
-
size_t getNMasks()¶
get the number of detector masks
-
void initBuffer(bool bufferAll = true)¶
Initialise the frame buffer.
-
void clearBuffer()¶
Clear the frame buffer.
-
inline bool allFramesBuffered() const¶
Check whether all frames are buffered.
Public Members
-
std::size_t datashape[3] = {0, 0, 0}¶
Data shape (columns, rows, frames)
-
std::size_t nFrames() const¶
-
class ImageGradient¶
- #include <ImageGradient.h>
Compute image gradients.
Compute image gradients via real space operations or convolution. Real space method is provided mainly as a sanity check for FFT method.
Public Functions
-
void reset()¶
Reset the stored matrices.
-
double pixel(int row, int col)¶
Get an element, respecting periodic boundary condition.
-
inline Eigen::MatrixXd *dx()¶
Get the gradient in the x direction.
-
inline Eigen::MatrixXd *dy()¶
Get the gradient in the y direction.
-
void compute(GradientKernel kernel, bool realspace)¶
Compute the gradient (but not the magnitude)
-
void computeRealSpace(GradientKernel kernel)¶
Compute gradient using real space method.
-
void computeFFT(GradientKernel kernel)¶
Compute gradient using FFT method.
-
void gradient(std::function<void(int, int)> kernel_operator)¶
Loop over the image to compute the gradient the gradient.
-
void centralDifference(int row, int col)¶
Compute gradient of a pixel with central differences.
-
void sobel(int row, int col)¶
Compute gradient of a pixel with Sobel operator.
-
void prewitt(int row, int col)¶
Compute gradient of a pixel with Prewitt operator.
-
void roberts(int row, int col)¶
Compute gradient of a pixel with Roberts operator.
-
Eigen::MatrixXd magnitude() const¶
Return the magnitude of the gradient.
-
void reset()¶
-
class SingleFrame : public ohkl::DataSet¶
- #include <SingleFrame.h>
Special case of DataSet that contains only one image.
Public Functions
-
void addFrame(const std::string &rawfilename, DataFormat format)¶
Add a raw file to be read as a single detector image frame. Only allow one frame to be added.
-
virtual Eigen::MatrixXi frame(const std::size_t idx) const override¶
Read a single frame.
-
virtual Eigen::MatrixXd transformedFrame(std::size_t idx) const override¶
Returns frame after transforming to account for detector gain and baseline.
-
virtual Eigen::MatrixXd gradientFrame(std::size_t idx, GradientKernel kernel, bool realspace = true) const override¶
Return per-pixel magnitude of gradient of a given frame.
-
std::size_t nFrames() const¶
Number of detector image frames acquired so far.
-
std::size_t nRows() const¶
The number of rows in each detector image.
-
std::size_t nCols() const¶
The number of columns in each detector image.
-
void addBoxMask(const AABB &aabb)¶
Add a rectangular mask (swig-compatible)
-
void addEllipseMask(const AABB &aabb)¶
Add an elliptical mask (swig-compatible)
-
void open()¶
Gets the file handle.
-
void close()¶
Close file and release handle.
-
ReciprocalVector computeQ(const DetectorEvent &ev) const¶
Returns the sample-space q vector corresponding to a detector event.
-
const IDataReader *reader() const¶
Returns the data reader used to set this dataset.
-
const Diffractometer *diffractometer() const¶
Returns the diffractometer associated to this dataset.
-
Diffractometer *diffractometer()¶
Returns the diffractometer associated to this dataset.
-
Detector &detector()¶
Returns the detector associated to this dataset.
-
const Detector &detector() const¶
Returns the detector associated to this dataset.
-
const ohkl::MetaData &metadata() const¶
Returns a const reference to the MetaData container.
-
ohkl::MetaData &metadata()¶
Returns a reference to the MetaData container.
-
void addDataFile(const std::string &filename, const DataFormat &fmt)¶
Add a data file for reading data. Reading frames will be done only upon request.
-
void setImageReaderParameters(const DataReaderParameters ¶ms)¶
Set the parameters for the image (raw or tiff) reader.
-
void addFrame(const std::string &filename, const DataFormat &format)¶
Add an image from a sequence.
-
void finishRead()¶
Finish reading procedure (must be called before using the data stored in the DataSet).
-
double wavelength() const¶
Query the wavelength stored in the metadata.
-
BitDepth bitDepth() const¶
Query image bit depth stored in metadata.
-
void setInstrumentStates(InstrumentStateSet *states)¶
Get the initial instrument states.
-
InstrumentStateList &instrumentStates()¶
Return instrument state list.
-
void adjustDirectBeam(double x_offset, double y_offset)¶
Adjust the direct beam (mainly for scripting purposes)
-
void initHistograms(std::size_t nbins)¶
Initialise intensity histograms.
-
void getFrameIntensityHistogram(std::size_t index)¶
Generate frame intensity histogram.
-
void clearHistograms()¶
Free histogram memeory.
-
inline size_t getNumberHistograms()¶
Get number of available histograms.
-
gsl_histogram *getHistogram(int index)¶
accessing created histograms
-
gsl_histogram *getTotalHistogram()¶
accessing Total histogram
-
bool hasMasks()¶
returns a booleans whether masks have been created or not
-
size_t getNMasks()¶
get the number of detector masks
-
void initBuffer(bool bufferAll = true)¶
Initialise the frame buffer.
-
void clearBuffer()¶
Clear the frame buffer.
-
inline bool allFramesBuffered() const¶
Check whether all frames are buffered.
Public Members
-
std::size_t datashape[3] = {0, 0, 0}¶
Data shape (columns, rows, frames)
-
void addFrame(const std::string &rawfilename, DataFormat format)¶
-
class Experiment¶
- #include <Experiment.h>
Top level core object wrapping all other core objects.
The Experiment object contains organises and controls almost all other core classes. It is a container for DataSets, PeakCollections and UnitCells. These objects are managed via their handlers, DataHandler, PeakHandler and CellHandler. It also calls the principal algorithms for finding peaks, integration, autoindexing, refining and merging peak collections.
Public Functions
-
const std::string &name() const¶
Get the name of the Experiment.
-
void setName(const std::string &name)¶
Set the name of the Experiment.
-
void readFromYaml(const std::string &filename)¶
Read saved experiment parameters from yaml.
-
void saveToYaml(const std::string &filename)¶
Save the parameters of the experiment to yaml.
-
Diffractometer *getDiffractometer()¶
Get a pointer to the diffractometer.
-
void setDiffractometer(const std::string &diffractometerName)¶
Set the named diffractometer from the relevant YAML instrument file.
-
inline DataReaderParameters *dataReaderParameters() const¶
Get a pointer to the data reader parameters.
-
const DataMap *getDataMap() const¶
Get the map of data sets from the handler.
-
sptrDataSet getData(const std::string &name) const¶
Get shared pointer to named data set.
-
DataList getAllData()¶
Get a vector containing shared pointers to all data sets.
-
int numData() const¶
Return number of data sets.
-
bool addData(sptrDataSet data, bool default_states = true)¶
Add a data set to the handler.
-
bool hasData(const std::string &name) const¶
Check whether the handler has named data set.
-
void removeData(const std::string &name)¶
Remove data set from handler.
-
bool addPeakCollection(const std::string &name, const PeakCollectionType type, std::vector<Peak3D*> peaks, sptrDataSet data, sptrUnitCell cell)¶
Create a new PeakCollection from a vector of peaks.
-
bool hasPeakCollection(const std::string &name)¶
Check if the handler has the named peak collection.
-
bool hasPeakCollectionType(PeakCollectionType t)¶
Check if handler has Peak Collections of a certain type.
-
bool hasIntegratedPeakCollection()¶
returns whether or not experiment has integrated peakcollections
-
PeakCollection *getPeakCollection(const std::string name)¶
Get a pointer to the named peak collection.
-
void removePeakCollection(const std::string &name)¶
Remove the named peak collection.
-
int numPeakCollections() const¶
Get the number of peak collections.
-
bool acceptFilter(std::string name, PeakCollection *collection, PeakCollectionType pct, sptrDataSet data)¶
Create a new peak collection from peaks caught by a filter.
-
bool clonePeakCollection(std::string name, std::string new_name)¶
Duplicate a peak collection (deep copy) for comparison after some process.
-
void addUnitCell(const std::string &name, const UnitCell &unit_cell)¶
Add a unit cell to the experiment.
-
void addUnitCell(const std::string &name, double a, double b, double c, double alpha, double beta, double gamma, sptrDataSet data)¶
Add a unit cell to the experiment via cell parameters (skip autoindexing step)
-
void addUnitCell(const std::string &name, double a, double b, double c, double alpha, double beta, double gamma, const std::string &space_group, sptrDataSet data)¶
Add a user-defined unit cell including space group.
-
bool hasUnitCell(const std::string &name) const¶
Returns true if the experiment has a data.
-
std::vector<std::string> getUnitCellNames() const¶
Get a list of loaded list names.
-
UnitCell *getUnitCell(const std::string &name) const¶
Return a pointer to the named unit cell.
-
sptrUnitCell getSptrUnitCell(const std::string &name) const¶
Return a pointer to the named unit cell.
-
sptrUnitCell getSptrUnitCell(const unsigned int id) const¶
Return a pointer to the numbered unit cell.
-
void removeUnitCell(const std::string &name)¶
Remove a unit cell from the experiment.
-
void swapUnitCells(const std::string &old_cell, const std::string &new_cell)¶
Swap two unit cells in the map contained by the handler.
-
int numUnitCells() const¶
Get the number of unit cells in the map.
-
bool checkAndAssignUnitCell(PeakCollection *peaks, double length_tol, double angle_tol, std::string name = kw_acceptedUnitcell)¶
Check solution against reference cell and accept if within tolerances.
-
void assignUnitCell(PeakCollection *peaks, std::string cellName = ohkl::kw_acceptedUnitcell)¶
Assign unit cell to a peak collection, compute Miller indices from q and cell.
-
void setReferenceCell(double a, double b, double c, double alpha, double beta, double gamma, sptrDataSet data)¶
Set the reference cell.
-
std::vector<std::string> getCompatibleSpaceGroups() const¶
Get space groups compatible with unit cell.
-
UnitCellHandler *getCellHandler() const¶
Get the cell handler.
-
void setLastUnitCellIndex(unsigned int index)¶
set last unit cell index in cell handler
-
std::vector<UnitCell*> getUnitCells()¶
Get a vector of unit cells in the experiment.
-
std::vector<sptrUnitCell> getSptrUnitCells(sptrDataSet data = nullptr)¶
Get a vector of unit cells in the experiment.
-
std::vector<PeakCollection*> getPeakCollections(sptrDataSet data = nullptr)¶
get a vector of pointers to peak collections
-
bool addShapeModel(const std::string &name, const ShapeModel &shapes)¶
Add a shape model.
-
bool addShapeModel(const std::string &name, std::unique_ptr<ShapeModel> &shapes)¶
Add a shape model from a unique_ptr.
-
bool hasShapeModel(const std::string &name) const¶
Returns true if the experiment has named shape model.
-
ShapeModel *getShapeModel(const std::string name)¶
Returns the named shape model.
-
void removeShapeModel(const std::string &name)¶
Remove a shape model from the experiment.
-
int numShapeModels() const¶
Get the number of shape models.
-
std::string generateShapeModelName()¶
Generate name for new shape model.
-
std::vector<ShapeModel*> getShapeModels(sptrDataSet data = nullptr)¶
Get a vector of pointers to shape models.
-
bool addInstrumentStateSet(sptrDataSet data)¶
Add a set of instrment states.
-
bool addInstrumentStateSet(sptrDataSet data, const InstrumentStateList &states, bool overwrite = true)¶
Add a set of instrment states.
-
bool addInstrumentStateSet(sptrDataSet data, std::unique_ptr<InstrumentStateSet> &states)¶
Add a set of instrment states.
-
InstrumentStateSet *getInstrumentStateSet(const sptrDataSet &data)¶
Returns the named InstrumentStateSet.
-
InstrumentStateSet *getInstrumentStateSet(const DataSet *data)¶
Returns the named InstrumentStateSet.
-
void removeInstrumentStateSet(const sptrDataSet &data)¶
Remove a set of instrument states from the experiment.
-
int numInstrumentStateSets() const¶
Return number of instrument state sets.
-
inline PeakFinder *peakFinder()¶
Return a pointer to the PeakFinder object.
-
inline PeakFinder2D *peakFinder2D()¶
Return a pointer to the OpenCV peak finder.
-
bool acceptFoundPeaks(const std::string &name)¶
Create a new peak collection from the peaks found by the peak finder.
-
bool acceptFoundPeaks(const std::string &name, const PeakCollection &found)¶
Create a new peak collection from a found collection.
-
inline PeakFilter *peakFilter()¶
Return a pointer to the PeakFilter object.
-
inline AutoIndexer *autoIndexer() const¶
Get a pointer to the AutoIndexer object.
-
const UnitCell *getAcceptedCell() const¶
Get a pointer to the accepted/assigned unit cell.
-
const UnitCell *getReferenceCell() const¶
Get a pointer to the reference unit cell.
-
inline ShapeModelBuilder *shapeModelBuilder() const¶
Get a pointer to the ShapeModelBuilder object.
-
inline PeakMerger *peakMerger() const¶
get a pointer to the PeakMerger object
-
void saveToFile(const std::string &path) const¶
Save the current experiment state to hdf5.
-
void loadFromFile(const std::string &path)¶
Load the current experiment state to hdf5.
-
Integrator *integrator()¶
Get a pointer to the integrator module.
-
std::string generatePeakCollectionName()¶
Generate automatic name for PeakCollection.
-
std::string generateUnitCellName()¶
Generate automatic name for UnitCell.
-
inline void setStrategy(bool flag)¶
Toggle strategy mode.
-
inline bool strategy() const¶
Get strategy mode flag.
-
const std::string &name() const¶
-
class ExperimentYAML¶
- #include <ExperimentYAML.h>
Read and write yml file containing experiment parameters.
Implements yaml parser for yml file containing experiment parameters (NOT) data or reduced data). Mainly for GUI purposes, but may have other uses.
Public Functions
-
void grabDataReaderParameters(DataReaderParameters *params) const¶
Get parameters from data reader.
-
void setDataReaderParameters(DataReaderParameters *params)¶
Set data reader parameters.
-
void grabIntegrationParameters(IntegrationParameters *params)¶
Get parameters from integrator.
-
void setIntegrationParameters(IntegrationParameters *params)¶
Set integrator parameters.
-
void grabPeakFinderParameters(PeakFinderParameters *params)¶
Get parameters from integrator.
-
void setPeakFinderParameters(PeakFinderParameters *params)¶
Set peak finder parameters.
-
void grabAutoindexerParameters(IndexerParameters *params)¶
Get parameters from indexer.
-
void setAutoindexerParameters(IndexerParameters *params)¶
Set indexer parameters.
-
void grabShapeParameters(ShapeModelParameters *params)¶
Get parameters from shape model.
-
void setShapeParameters(ShapeModelParameters *params)¶
Set shape model parameters.
-
void grabPredictorParameters(PredictionParameters *params)¶
Get parameters from predictor.
-
void setPredictorParameters(PredictionParameters *params)¶
Set shape model parameters.
-
void grabMergeParameters(MergeParameters *params)¶
Get parameters from merger.
-
void setMergeParameters(MergeParameters *params)¶
Set merger parameters.
-
void writeFile(const std::string &filename)¶
Write the yml file.
-
void grabDataReaderParameters(DataReaderParameters *params) const¶
-
class Integrator¶
- #include <Integrator.h>
Integrate predicted peaks.
Implements naive pixel sum integration and various profile fitting integration methods for computing intensities, sigmas and strengths of reflections.
Construct a shape collection of integrated strong peaks, and compute the mean covariance of shapes in the vicinity of a predicted peaks to generate an integration region. This is integrated via profile integration.
Public Functions
-
ohkl::IIntegrator *getIntegrator(const IntegratorType integrator_type) const¶
Get an integrator from the map.
-
DataHandler *getDataHandler()¶
Return a pointer to the data handler.
-
void integratePeaks(sptrDataSet data, PeakCollection *peaks, IntegrationParameters *params, ShapeModel *shapes, bool parallel = true)¶
Integrate a peak collection.
-
void integrateFoundPeaks(PeakFinder *peak_finder, bool parallel = true)¶
Integrate peaks found by _peak_finder.
Set the parameters.
-
IntegrationParameters *parameters()¶
Get the parameters.
-
void setHandler(sptrProgressHandler handler)¶
Set a progress handler.
-
unsigned int numberOfPeaks()¶
Get the number of valid peaks;.
-
unsigned int numberOfValidPeaks()¶
Get the total number of peaks.
-
ohkl::IIntegrator *getIntegrator(const IntegratorType integrator_type) const¶
-
struct PeakFinderParameters¶
- #include <PeakFinder.h>
Parameters for finding peaks.
Public Members
-
int minimum_size = 30¶
Minimum number of pixels in a blob.
-
int maximum_size = 10000¶
Maximum number of pixels in a blob.
-
double peak_end = 1.0¶
Peak scaling factor (sigmas)
-
int maximum_frames = 10¶
Maximum number of frames a blob can span.
-
int first_frame = -1¶
First frame in range for peak finding.
-
int last_frame = -1¶
Last frame in range for peak findinng.
-
double threshold = 80.0¶
Blobs containing fewer counts than this are discarded.
-
int minimum_size = 30¶
-
class PeakFinder¶
- #include <PeakFinder.h>
Perform image recognition on detector images to find peaks.
Given a “DataSet” (set of detector images from one experimental run), PeakFinder will locate all detector spots per-frame, and merge sets of contiguous partial reflections into a blog in 3D space (a “Blob3D”). A blob will become a Bragg reflexis for which an intensity and sigma can be calculated.
Public Functions
-
PeakFinder()¶
Default constructor.
-
void find(const sptrDataSet data)¶
Iterate through a vector of DataSets, finding the peaks for each.
-
std::vector<Peak3D*> currentPeaks()¶
Return a vector of peaks found by PeakFinder.
-
PeakCollection *getPeakCollection()¶
Return a pointer to the PeakCollection owned by PeakFinder.
Generate a PeakCollection from a vector of found peaks.
-
inline sptrDataSet currentData()¶
Return the DataList (vector of pointers to DataSets)
-
void setHandler(const sptrProgressHandler &handler)¶
Set the progress handler.
-
PeakFinderParameters *parameters()¶
Get the parameters for peak finding.
-
unsigned int numberFound()¶
Get the number of peaks found.
-
inline bool foundPeaks() const¶
Have we found peaks?
-
inline bool isIntegrated()¶
Return integrated state.
-
inline bool hasBkgGradient()¶
Return background gradient status.
-
inline void setIntegrated(bool b)¶
set integration state
-
inline void setBkgGradient(bool b)¶
set background gradient status
-
void setConvolver(std::unique_ptr<Convolver> convolver)¶
Set the convolver flavour for peak/background convolution.
-
void setConvolver(const std::string &convolver, const std::map<std::string, double> ¶meters)¶
Set the convolver flavour for peak/background convolution.
-
void setConvolver(const Convolver &convolver)¶
Set the convolver flavour for peak/background convolution.
-
inline ohkl::Convolver *convolver() const¶
Get the convolver.
-
PeakFinder()¶
-
struct PeakFinder2DParameters : public cv::SimpleBlobDetector::Params¶
- #include <PeakFinder2D.h>
Perform image recognition on detector images to find peaks in 2D.
PLACEHOLDER
-
class PeakFinder2D¶
- #include <PeakFinder2D.h>
Public Functions
-
void find(std::size_t image_idx)¶
Find blobs on given image.
-
void findAll()¶
Find blobs for all images in the data set.
-
void setHandler(const sptrProgressHandler &handler)¶
Set the progress handler.
-
inline PeakFinder2DParameters *parameters()¶
Get the parameters.
-
inline KeyPointCollection *keypoints()¶
Get object storing per-frame keypoints.
-
void setConvolver(const ConvolutionKernelType &kernel)¶
Set the convolver.
-
inline Convolver *convolver() const¶
Get the convolver.
-
void find(std::size_t image_idx)¶
-
class ShapeModelBuilder¶
- #include <ShapeModelBuilder.h>
Build a shape model.
-
class Diffractometer¶
- #include <Diffractometer.h>
#include “core/loader/IDataReader.h” diffractometer setup, consisting of a sample, source, and detector.
Public Functions
-
const std::string &name() const¶
Returns the name of this diffractometer.
-
void setName(const std::string &name)¶
Sets the name of the diffractometer.
-
Detector *detector()¶
Returns a pointer to the detector of this diffractometer.
-
const Detector *detector() const¶
Returns const pointer to the detector of this diffractometer.
-
InstrumentState instrumentState(const std::size_t frame_idx)¶
Returns the instrument state.
-
void setDetector(std::unique_ptr<Detector> detector)¶
Sets the detector of this diffractometer.
-
Sample &sample()¶
Returns the non-const reference to the sample of this diffractometer.
-
const Sample &sample() const¶
Returns the const reference to the sample of this diffractometer.
-
void setSample(const Sample &sample)¶
Sets the sample of this diffractometer.
-
Source &source()¶
Returns the non-const reference to the source of this diffractometer.
-
const Source &source() const¶
Returns the const reference to the source of this diffractometer.
-
void setSource(const Source &source)¶
Sets the source of this diffractometer.
-
void addSampleAngles(std::size_t frame_idx, const DataReaderParameters ¶ms)¶
Add angles for a sample state.
-
void addSampleAngles(const std::vector<double> &angles)¶
Add angles for a given sample state.
-
void addDetectorAngles(const DataReaderParameters ¶ms)¶
Add angles for a detector state.
-
void addDetectorAngles(const std::vector<double> &angles)¶
Add angles for a detector state.
-
inline const std::vector<std::vector<double>> &sampleAngles() const¶
Get the sample angles.
-
inline const std::vector<std::vector<double>> &detectorAngles() const¶
Get the detector angles.
-
void setSampleAngles(const RowMatrixXd &mat, std::size_t nframes)¶
Set the sample angles.
-
void setDetectorAngles(const RowMatrixXd &mat, std::size_t nframes)¶
Set the detector angles.
-
const std::string &name() const¶
-
class InstrumentState¶
- #include <InstrumentState.h>
Class storing the state of the experiment at a given point in time.
State refers to any parameters which might change during the experiment: sample orientation, sample position, etc. States are initially loaded as metadata but can also be refined as part of the data treatment. Note that InstrumentState can be invalid (_valid = false) if interpolation fails
Subclassed by ohkl::InterpolatedState
Public Functions
-
InstrumentState(Diffractometer *diffractometer = nullptr)¶
default value needed for SWIG (note: nullptr does not work)
-
ReciprocalVector kfLab(const DirectVector &detector_position) const¶
Takes a direct vector in detector coordinates and returns kf in lab coordinates.
-
ReciprocalVector ki() const¶
Returns source wavevector k_i.
-
void adjustKi(const DirectVector &detector_position)¶
Adjust ki by shifting the direct beam on the detector.
-
ReciprocalVector sampleQ(const DirectVector &detector_position) const¶
Takes direct vector in detector coordinates and returns q in sample coordinates.
-
double gamma(const DirectVector &detector_position) const¶
Returns the gamma angle associated to the given lab space position.
-
double nu(const DirectVector &detector_position) const¶
Returns the nu angle associated to the given lab space position.
-
double twoTheta(const DirectVector &detector_position) const¶
Returns the 2*theta angle associated to the given lab space position.
-
Eigen::Matrix3d jacobianK(double px, double py) const¶
Compute the jacobian of the transformation (x,y) -> k_lab.
-
void setDiffractometer(Diffractometer *diffractometer)¶
Set the diffractometer.
-
Diffractometer *diffractometer()¶
Returns a pointer to the diffractometer of the state.
-
const Diffractometer *diffractometer() const¶
Returns a const pointer to the diffractometer of the state.
-
bool isValid() const¶
Returns whether the InstrumentState is valid.
-
Eigen::Matrix3d detectorOrientationMatrix() const¶
compute the sample orientation from fixed orientation and offset
-
Eigen::Matrix3d sampleOrientationMatrix() const¶
compute the sample orientation from fixed orientation and offset
-
std::string toString() const¶
Return state as a string.
Public Members
-
Eigen::Matrix3d detectorOrientation¶
Detector orientation as read from metadata.
-
Eigen::Quaterniond sampleOrientation¶
Sample orientation as read from metadata.
-
Eigen::Quaterniond sampleOrientationOffset¶
Offset to sample orientation, used for parameter refinement.
-
Eigen::Vector3d samplePosition¶
Sample position.
-
Eigen::Vector3d detectorPositionOffset¶
Detector offset.
-
Eigen::RowVector3d ni¶
Incoming beam direction.
-
double wavelength¶
Incoming beam wavelength.
-
bool refined¶
True if this state has been refined.
Public Static Functions
-
static InstrumentState state(Diffractometer *const diffractometer, const std::size_t frame_idx)¶
Returns the instrument state as read from the metadata.
-
InstrumentState(Diffractometer *diffractometer = nullptr)¶
-
class InstrumentStateSet¶
- #include <InstrumentStateSet.h>
Container for mutable instrument states (sample positions/orientation, detector offset and incident wavevector) that are modified during refinement.
Public Functions
-
inline std::string name() const¶
Get the name of the InstrumentStateSet.
-
inline unsigned int id() const¶
Get the integer id.
-
inline void setId(unsigned int id)¶
Set the integer id.
-
inline void reset()¶
Delete the instrument states.
-
inline void setInstrumentStates(const InstrumentStateList &states)¶
Set the states from a vector of InstrumentStates.
-
inline InstrumentStateList getInstrumentStateList()¶
Get the vector of InstrumentStates.
-
const InstrumentState *state(std::size_t frame)¶
Return a single InstrumentState by frame index.
-
void applyBeamOffset()¶
Apply the direct beam offset specified by the instrument (yml2c) file.
-
inline std::string name() const¶
-
class InterpolatedState : public ohkl::InstrumentState¶
- #include <InterpolatedState.h>
Stores an interpolation between two instrument states.
In addition to the data fields stored in an InstrumentState, this class also stores information about the (angular) velocity and step size, which is needed to compute analytic derivatives of various functions.
Public Functions
-
InterpolatedState(Diffractometer *diffractometer = nullptr)¶
Default value needed for SWIG (note: nullptr does not work)
-
InterpolatedState(const InstrumentState &s1, const InstrumentState &s2, double t, bool valid = true)¶
Construct by interpolation. The paramter t should be between 0 and 1.
-
Eigen::Matrix3d jacobianQ(double px, double py) const¶
Compute the jacobian of the transformation (x,y,frame) -> q_sample.
-
double lorentzFactor(double px, double py) const¶
Compute the Lorentz factor at the detector coordinates (px, py).
-
ReciprocalVector kfLab(const DirectVector &detector_position) const¶
Takes a direct vector in detector coordinates and returns kf in lab coordinates.
-
ReciprocalVector ki() const¶
Returns source wavevector k_i.
-
void adjustKi(const DirectVector &detector_position)¶
Adjust ki by shifting the direct beam on the detector.
-
ReciprocalVector sampleQ(const DirectVector &detector_position) const¶
Takes direct vector in detector coordinates and returns q in sample coordinates.
-
double gamma(const DirectVector &detector_position) const¶
Returns the gamma angle associated to the given lab space position.
-
double nu(const DirectVector &detector_position) const¶
Returns the nu angle associated to the given lab space position.
-
double twoTheta(const DirectVector &detector_position) const¶
Returns the 2*theta angle associated to the given lab space position.
-
Eigen::Matrix3d jacobianK(double px, double py) const¶
Compute the jacobian of the transformation (x,y) -> k_lab.
-
void setDiffractometer(Diffractometer *diffractometer)¶
Set the diffractometer.
-
Diffractometer *diffractometer()¶
Returns a pointer to the diffractometer of the state.
-
const Diffractometer *diffractometer() const¶
Returns a const pointer to the diffractometer of the state.
-
bool isValid() const¶
Returns whether the InstrumentState is valid.
-
Eigen::Matrix3d detectorOrientationMatrix() const¶
compute the sample orientation from fixed orientation and offset
-
Eigen::Matrix3d sampleOrientationMatrix() const¶
compute the sample orientation from fixed orientation and offset
-
std::string toString() const¶
Return state as a string.
Public Members
-
Eigen::Vector3d axis¶
The axis of crystal rotation, in sample space.
-
double stepSize¶
Step size between the two underlying InstrumentStates, in radians.
-
Eigen::Matrix3d detectorOrientation¶
Detector orientation as read from metadata.
-
Eigen::Quaterniond sampleOrientation¶
Sample orientation as read from metadata.
-
Eigen::Quaterniond sampleOrientationOffset¶
Offset to sample orientation, used for parameter refinement.
-
Eigen::Vector3d samplePosition¶
Sample position.
-
Eigen::Vector3d detectorPositionOffset¶
Detector offset.
-
Eigen::RowVector3d ni¶
Incoming beam direction.
-
double wavelength¶
Incoming beam wavelength.
-
bool refined¶
True if this state has been refined.
Public Static Functions
-
static InterpolatedState interpolate(const InstrumentStateList &states, const double frame_idx)¶
Interpolates the state for a given frame index (not necessarily integer)
-
static InstrumentState state(Diffractometer *const diffractometer, const std::size_t frame_idx)¶
Returns the instrument state as read from the metadata.
-
InterpolatedState(Diffractometer *diffractometer = nullptr)¶
-
class Blob3D¶
- #include <Blob3D.h>
A region of interest in a 3D image.
A Blob is constructed by adding points in the image with coordinates x,y,z and an associated mass that represents any scalar field such as intensity. Blob objects records the total mass, the mass-weighted first and second moments as new points are added to the blob. Knowledge about individual points is lost, i.e Blob can only increase in size.Blobs can be merged and maintain zero, first and second momentum. Blob3D can be transformed into an Ellipsoid, by diagonalizing the variance tensor.
Public Functions
-
Blob3D()¶
Initialize an empty blob.
-
Blob3D(double x, double y, double z, double m)¶
Initialize a blob with a point of mass m at x,y,z.
-
void addPoint(double x, double y, double z, double m)¶
Add point to the Blob.
-
double getMass() const¶
Returns the total mass.
-
int getComponents() const¶
Returns the number of points.
-
double getMinimumMass() const¶
Returns the minimumMass.
-
double getMaximumMass() const¶
Returns the minimumMass.
-
Eigen::Vector3d center() const¶
Returns the center of Mass.
-
bool toEllipsoid(double scale, Eigen::Vector3d ¢er, Eigen::Vector3d &eigenvalues, Eigen::Matrix3d &eigenvectors) const¶
Gets the ellipsoid parameters.
-
void printSelf(std::ostream &os) const¶
Print to ostream.
-
Eigen::Matrix3d covariance() const¶
Gets covariance matrix of the blob.
-
bool isValid() const¶
Whether the blob is valid.
-
Blob3D()¶
-
class GaussianIntegrator : public ohkl::IIntegrator¶
- #include <GaussianIntegrator.h>
Compute integrated intensity by fitting to an analytic 3D Gaussian.
Public Functions
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Integrate a vector of peaks.
-
void parallelIntegrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Thread-parallel integrate.
-
void setHandler(sptrProgressHandler handler)¶
Set the progress handler.
-
void setParameters(const IntegrationParameters ¶ms)¶
Assign a parameter set to the integrator.
-
inline void setParallel(bool parallel)¶
Toggle parallel integration.
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
-
struct ComputeResult¶
- #include <IIntegrator.h>
Return value of all integrator compute methods. Contains any mutable data from integration.
Public Functions
-
ComputeResult()¶
Construct an invalid result by default.
Public Members
-
RejectionFlag integration_flag¶
Rejection flag for integration.
-
Intensity sum_intensity¶
Pixel sum integrated intensity after background subtraction.
-
Intensity profile_intensity¶
Profile integrated intensity after background correction.
-
Intensity sum_background¶
Mean local pixel sum background of peak. The uncertainty is the uncertainty of the estimate of the background.
-
Intensity profile_background¶
Profile integrated background of peak. The uncertainty is the uncertainty of the estimate of the background.
-
Intensity bkg_gradient¶
Mean gradient of background (Gaussian statistics)
-
std::vector<Intensity> rocking_curve¶
The rocking curve of the peak.
-
IntegratorType integrator_type¶
Type of integrator used.
-
std::optional<Ellipsoid> shape¶
Shape of peak when centre/covariance are changed during integration.
-
ComputeResult()¶
-
struct IntegrationParameters¶
- #include <IIntegrator.h>
Structure containing parameters for all integrators.
Subclassed by ohkl::PredictionParameters, ohkl::ShapeModelParameters
Public Members
-
double peak_end = 3.0¶
End of peak region for RegionType::VariableEllipsoid (sigmas)
-
double bkg_begin = 3.0¶
Beginning of background region (sigmas)
-
double bkg_end = 6.0¶
End of background region (sigmas)
-
double fixed_peak_end = 8¶
End of peak region for RegionType::FixedEllipsoid (pixels)
-
double fixed_bkg_begin = 1.0¶
Beginning of background region (factor of fixed_peak_end)
-
double fixed_bkg_end = 2.0¶
End of background region (factor of fixed_peak_end)
-
double max_counts = 50000.0¶
Maximum per-pixel count.
-
double max_strength = 1¶
Maximum strength for profile integration.
-
double max_d = 2.5¶
Maximum d (minimum resolution) for profile integration.
-
int max_width = 10¶
Maximum number of images a peak/shape can span.
-
bool fit_center = true¶
Whether to update the peak centre after integration.
-
bool fit_cov = true¶
Whether to update the peak covariance after integration.
-
bool discard_saturated = false¶
Whether to discard peaks with saturated pixels.
-
IntegratorType integrator_type = IntegratorType::PixelSum¶
Type of integrator.
-
bool use_gradient = false¶
Use gradient to discriminate heterogeneous background regions.
-
GradientKernel gradient_type = GradientKernel::Sobel¶
Kernel to use for gradient convolution.
-
bool fft_gradient = false¶
Whether to use FFT or real space gradient computation.
-
RegionType region_type = RegionType::VariableEllipsoid¶
Whether to use fixed or sigma-dependent integration regions.
-
bool skip_masked = true¶
Whether to skip peaks that intersect masks.
-
bool remove_overlaps = false¶
Whether to remove peaks with overlapping peak integration regions.
-
bool use_max_strength = false¶
Skip profile integration of peaks with strength above threshold.
-
bool use_max_d = false¶
Skip profile integration of peaks with resolution below threshold.
-
bool use_max_width = false¶
Skip integration of peaks spanning too many images.
-
double peak_end = 3.0¶
-
class IIntegrator¶
- #include <IIntegrator.h>
Base class for integrators. Handles per-frame integration of a peak.
All integrators inherit from this class.
Subclassed by ohkl::GaussianIntegrator, ohkl::PixelSumIntegrator, ohkl::Profile1DIntegrator, ohkl::Profile3DIntegrator
Public Functions
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Integrate a vector of peaks.
-
void parallelIntegrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Thread-parallel integrate.
-
void setHandler(sptrProgressHandler handler)¶
Set the progress handler.
-
void setParameters(const IntegrationParameters ¶ms)¶
Assign a parameter set to the integrator.
-
inline void setParallel(bool parallel)¶
Toggle parallel integration.
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
-
class ISigmaIntegrator : public ohkl::PixelSumIntegrator¶
- #include <ISigmaIntegrator.h>
Compute the integrated intensity via the I/sigma method used in RETREAT.
Described in E. Prince et al., J. Appl. Cryst., 30:133, 1997. doi:10.1107/S0021889896012824.
We take the numeric minimum,
\[ \textrm{argmin}\left( \frac{\sigma^2(p)}{p^2} + \frac{\sigma^2(I)}{I^2} \right), \]Which determines the end of the peak region, and use this to do profile fitting integration.Public Functions
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Integrate a vector of peaks.
-
void parallelIntegrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Thread-parallel integrate.
-
void setHandler(sptrProgressHandler handler)¶
Set the progress handler.
-
void setParameters(const IntegrationParameters ¶ms)¶
Assign a parameter set to the integrator.
-
inline void setParallel(bool parallel)¶
Toggle parallel integration.
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
-
class PixelSumIntegrator : public ohkl::IIntegrator¶
- #include <PixelSumIntegrator.h>
Peak integration using naive background estimation and subtraction.
Subclassed by ohkl::ISigmaIntegrator, ohkl::ShapeIntegrator
Public Functions
-
PixelSumIntegrator()¶
Construct the pixel sum integrator
- Parameters:
fit_center – update the peak center as part of integration
fit_covariance – update the peak shape covariance matrix as part of integration
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Integrate a vector of peaks.
-
void parallelIntegrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Thread-parallel integrate.
-
void setHandler(sptrProgressHandler handler)¶
Set the progress handler.
-
void setParameters(const IntegrationParameters ¶ms)¶
Assign a parameter set to the integrator.
-
inline void setParallel(bool parallel)¶
Toggle parallel integration.
-
PixelSumIntegrator()¶
-
class Profile1DIntegrator : public ohkl::IIntegrator¶
- #include <Profile1DIntegrator.h>
Peak integrator using 1D profile fitting.
Described in:
W. Kabsch, J. Appl. Crystallography, 21:916, 1988. doi:10.1107/S0021889888007903
W. Kabsch, Acta Crystallographica D, 66:133, 2010. doi:10.1107/s0907444909047374
The implementation is identical to the Profile3DIntegrator, except the susbscripts \(i\) refer to elements of a 1D profile instead of a 3D profile.
Public Functions
-
Profile1DIntegrator()¶
Construct integrator.
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Integrate a vector of peaks.
-
void parallelIntegrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Thread-parallel integrate.
-
void setHandler(sptrProgressHandler handler)¶
Set the progress handler.
-
void setParameters(const IntegrationParameters ¶ms)¶
Assign a parameter set to the integrator.
-
inline void setParallel(bool parallel)¶
Toggle parallel integration.
-
class Profile3DIntegrator : public ohkl::IIntegrator¶
- #include <Profile3DIntegrator.h>
Integrate a peak using 3D profile fitting.
Described in:
W. Kabsch, J. Appl. Crystallography, 21:916, 1988. doi:10.1107/S0021889888007903
W. Kabsch, Acta Crystallographica D, 66:133, 2010. doi:10.1107/s0907444909047374
Here, we minimise the chi-squared loss subject to the normalisation \(\sum_i p_i^2 = 1\), Resulting in the 2x2 linear system:
\[\begin{split} \begin{bmatrix} \sum 1/v_i^2 && \sum p_i/v_i^2 \\ p_i/v_i^2 && p_i^2/v_i^2 \end{bmatrix} \begin{bmatrix} B \\ I \end{bmatrix} = \begin{bmatrix} \sum c_i/v_i^2 \\ \sum c_ip_i/v_i^2 \end{bmatrix} \end{split}\]Where \(c_i\), \(b_i\), \(v_i\) and \(p_i\) are the counts, background, variance and profile respectively, and I and B are the computed intensity and background respectively. These are solved via the following procedure:
Set \(v_i = b_i\) as an initial guess
Solve the linear equations by matrix inversion
Compute updated \(v_i = b_i - I p_i\)
Repeat from 2 until convergence
Public Functions
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Integrate a vector of peaks.
-
void parallelIntegrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Thread-parallel integrate.
-
void setHandler(sptrProgressHandler handler)¶
Set the progress handler.
-
void setParameters(const IntegrationParameters ¶ms)¶
Assign a parameter set to the integrator.
-
inline void setParallel(bool parallel)¶
Toggle parallel integration.
-
class ShapeIntegrator : public ohkl::PixelSumIntegrator¶
- #include <ShapeIntegrator.h>
Integrate a peak to generate a profile for ShapeModel.
Public Functions
-
ShapeIntegrator()¶
Construct the integrator with the given shape collection, bounding box, and box shape.
-
void initialise(const AABB &aabb, ShapeModelParameters *params)¶
Set shape model parameters.
-
void integrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Integrate a vector of peaks.
-
void parallelIntegrate(std::vector<Peak3D*> peaks, ShapeModel *shape_model, sptrDataSet data)¶
Thread-parallel integrate.
-
void setHandler(sptrProgressHandler handler)¶
Set the progress handler.
-
void setParameters(const IntegrationParameters ¶ms)¶
Assign a parameter set to the integrator.
-
inline void setParallel(bool parallel)¶
Toggle parallel integration.
-
ShapeIntegrator()¶
-
class Peak3D¶
- #include <Peak3D.h>
A peak object storing real-space information on the peak.
This object stores the real-space shape as an Ellipsoid object, including center of the detector spot in 3D (x, y, frame). The frame does not have to be an integer value, as the reflexis center is generally between two frames.
Also implements integration and reciprocal space transformation, as well as storing various metadata.
Public Functions
-
Peak3D(sptrDataSet dataSet)¶
Create peak belonging to data without setting a position, shape, or intensity.
-
Peak3D(sptrDataSet dataSet, const Ellipsoid &shape)¶
Create peak belonging to data with given shape.
-
Peak3D(sptrDataSet dataSet, const MillerIndex &hkl)¶
Create peak belonging to data with given shape.
Creat the peak from another peak.
-
void setShape(const Ellipsoid &shape)¶
Comparison operator used to sort peaks.
Sets the Peak region. Peak shaped is owned after setting
-
void setMillerIndices()¶
Set the Miller indices.
-
inline void setMillerIndices(const MillerIndex &hkl)¶
Set the Miller indices from the given triple. NB. ONLY TO BE USED BY EXPERIMENTIMPORTER.
-
const std::vector<Intensity> &rockingCurve() const¶
Gets the projection of total data in the bounding box.
-
Ellipsoid qShape() const¶
Compute the shape in q-space. May throw if there is no valid q-space ellipsoid.
This method computes an ellipsoid in q-space which is approximately the transformation from detector space to q-space of its shape ellipsoid (which is computed during blob search).
Suppose that the detector-space ellipsoid is given by the equation (x-x0).dot(A*(x-x0)) <= 1. Then if q = q0 + J(x-x0), then the corresponding ellipsoid.
This method can throw if there is no valid q-shape corresponding to the detector space shape.
-
const Ellipsoid &shape() const¶
Return the shape of the peak as an ellipsoid in detector coordinates.
-
Intensity correctedIntensity(const Intensity &intensity) const¶
Return intensity, after scaling, transmission, and Lorentz factor corrections.
-
Intensity correctedSumIntensity() const¶
Return corrected pixel sum intensity.
-
Intensity correctedProfileIntensity() const¶
Return corrected profile intensity.
-
inline Intensity sumIntensity() const¶
Return the pixel sum intensity of the peak.
-
inline Intensity profileIntensity() const¶
Return the profile integrated intensity of the peak.
-
inline Intensity sumBackground() const¶
Return pixel sum background of the peak.
-
inline Intensity profileBackground() const¶
Return the profile background.
-
inline Intensity meanBkgGradient() const¶
Return the mean background gradient of the peak.
-
inline void setSumIntensity(const Intensity &i)¶
Set the pixel sum intensity with no corrections.
-
inline void setProfileIntensity(const Intensity &i)¶
Set the profile integrated intensity with no corrections.
-
double peakEnd() const¶
Return shape scale used to define peak region.
-
double bkgBegin() const¶
Return shape scale used to define beginning of background region.
-
double bkgEnd() const¶
Return shape scale used to define end of background region.
-
double scale() const¶
Returns the scaling factor.
-
void setScale(double factor)¶
Set the scaling factor.
-
void reject(RejectionFlag flag)¶
Reject a peak.
-
bool enabled() const¶
Is the peak enabled (selected and not masked)?
-
void setTransmission(double transmission)¶
Set the transmission factor.
-
double transmission() const¶
Return the transmission factor.
-
void setUnitCell(const sptrUnitCell &uc)¶
Assign a unit cell to the peak.
-
const UnitCell *unitCell() const¶
Returns the unit cell.
-
void caughtYou(bool caught)¶
Designates the peak as “caught” by a filter.
-
void rejectYou(bool reject)¶
Designates the peak as “rejected” by a filter.
-
bool caughtByFilter() const¶
Has the peak been caught by a filter?
-
bool rejectedByFilter() const¶
Has the peak been rejected by a filter?
-
void setManually(const Intensity &sumInt, const Intensity &profInt, double peakEnd, double bkgBegin, double bkgEnd, int region_type, double scale, double transmission, const Intensity &sumBkg, const Intensity &profBkg, int rejection_flag, int sum_integration_flag, int profile_integration_flag, Intensity sumBkgGrad = {})¶
Manually set the integration parameters for this peak.
-
void updateIntegration(const ComputeResult &result, double peakEnd, double bkgBegin, double bkgEnd, const RegionType ®ionType)¶
Update the integration parameters for this peak.
-
ReciprocalVector q() const¶
Return the q vector of the peak, transformed into sample coordinates.
-
double d() const¶
Resolution (A) of this peak.
-
ReciprocalVector q(const InstrumentState &state) const¶
Return q vector in cases where we do not want to interpolate the InstrumentState (e.g. when indexing a single frame)
-
inline sptrDataSet dataSet() const¶
Return the data set to which this peak belongs.
-
const MillerIndex &hkl() const¶
Get the Miller indices for this peak.
-
inline double getPeakEnd()¶
Return the peak scale.
-
inline double getBkgBegin()¶
Return the beginniing of the background region (in peak scales)
-
inline double getBkgEnd()¶
Return the end of the background region (in peak scales)
-
inline RegionType regionType()¶
Return the integration region type.
-
void setRejectionFlag(const RejectionFlag &flag, bool overwrite = false)¶
Set the reason for this peak being disabled.
-
void setIntegrationFlag(const RejectionFlag &flag, const IntegratorType &integrator, bool overwrite = false)¶
Set the reason for rejection during integration.
-
void setMasked()¶
Label the peak as masked.
-
inline RejectionFlag getRejectionFlag() const¶
Return the rejection flag only.
-
inline RejectionFlag getSumIntegrationFlag() const¶
Return the sum integration flag only.
-
inline RejectionFlag getProfileIntegrationFlag() const¶
Return the profile integration flag only.
-
RejectionFlag rejectionFlag() const¶
Return the highest level rejection flag (pre-integration < sum < profile)
-
RejectionFlag sumRejectionFlag() const¶
Return the sum integration flag, or rejection flag if the former is not set.
-
RejectionFlag profileRejectionFlag() const¶
Return the sum integration flag, or rejection flag if the former is not set.
-
bool isRejectedFor(const RejectionFlag &flag) const¶
Check if any rejection flag matches argument.
-
std::string rejectionString() const¶
Return a string explaining the rejection.
-
std::string toString() const¶
Return a string representation of the peak.
-
void resetIntegration(IntegratorType integrator_type)¶
Reset integration state.
Public Static Functions
-
static const std::map<RejectionFlag, std::string> &rejectionMap()¶
Return the map of RejectionFlag definitions.
-
Peak3D(sptrDataSet dataSet)¶
-
class KeyPointCollection¶
- #include <KeyPointCollection.h>
Store a collection of OpenCV keypoints.
Used by PeakFinder2D to store blob centres, no shape information
-
class PeakCollection¶
- #include <PeakCollection.h>
Store a collection of peak shapes, for peak prediction and integration.
Container for a collection of peaks and its metadata
Public Functions
-
PeakCollection()¶
Default contructor.
-
inline unsigned int id()¶
Get integer id.
-
void setId(unsigned int id)¶
Set integer id.
-
void setName(const std::string &name)¶
Sets name of the PeakCollection.
-
std::string name() const¶
Returns the name of the PeakCollection.
-
inline sptrUnitCell unitCell() const¶
Get the unit cell from autoindexing or used to predict these peaks.
-
inline std::string description() const¶
Returns description of the PeakCollection.
-
inline void setDescription(std::string str)¶
set description
Populate the PeakCollection with a vector of shared pointers to peaks.
-
void populate(const std::vector<ohkl::Peak3D*> peak_list)¶
Populate the PeakCollection with a vector of raw pointers to peaks.
-
void push_back(const ohkl::Peak3D &peak)¶
Append one peak to the PeakCollection.
Append one peak to the PeakCollection.
-
void populateFromFiltered(PeakCollection *collection)¶
Populate from another collection, taking only peaks caughtByFilter.
-
void reset()¶
Remove all peaks from the PeakCollection.
-
std::vector<ohkl::Peak3D*> getFilteredPeakList() const¶
Return a std::vector of pointers to peaks caughtByFilter.
-
inline ohkl::PeakCollectionType type() const¶
Return the PeakCollectionType of the PeakCollection (FOUND, PREDICTED, etc.)
-
inline void setType(PeakCollectionType type)¶
Set the PeakCollectionType of the PeakCollection (FOUND, PREDICTED, etc.)
-
void computeSigmas()¶
compute beam divergence and mosaicity sigmas
-
double sigmaD() const¶
Return beam divergence sigma.
-
double sigmaM() const¶
Return mosaicity sigma.
-
void setMillerIndices() const¶
Set Miller indices of peaks for those with an assigned unit cell.
-
inline int numberOfPeaks() const¶
Get the number of peaks int he PeakCollection.
-
int numberOfValid() const¶
Return the number of valid peaks.
-
int numberOfInvalid() const¶
Return the number of invalid peaks.
-
MetaData &metadata()¶
Return a fresh generated pointer to metadata.
-
int countEnabled() const¶
Count enabled peaks.
-
inline bool isIndexed() const¶
Whether PeakCollection has been indexed or not.
-
inline bool isIntegrated() const¶
Whether PeakCollection has been integrated or not.
-
inline bool hasBkgGradient() const¶
Whether PeakCollection has background gradient.
-
inline void setIndexed(bool value)¶
Set Indexed flag.
-
inline void setIntegrated(bool value)¶
Count Integrated flag.
-
inline void setBkgGradient(bool value)¶
Background gradient flag.
-
void resetIntegration(IntegratorType integrator_type)¶
Reset integration status.
-
void resetRejectionFlags()¶
Reset peak rejection flags.
-
void resetIntegrationFlags(IntegratorType integrator)¶
Reset peak rejection status to pre-integration values.
-
PeakCollection()¶
-
struct PeakFilterFlags¶
- #include <PeakFilter.h>
Enable/disable the various types of filters.
Public Members
-
bool enabled¶
filter by enabled
-
bool masked¶
filter by mask
-
bool indexed¶
filter by indexed peaks
-
bool strength¶
filter by strength (intensity/sigma)
-
bool d_range¶
filter by detector range
-
bool extinct¶
filter by extinction (allowed by unit cell)
-
bool frames¶
catch peaks in a specifed frame range
-
bool rejection_flag¶
catch peaks with a specific rejection flag
-
bool intensity¶
filter by intensity
-
bool sigma¶
filter by sigma
-
bool gradient¶
filter by background gradient
-
bool gradient_sigma¶
filter by background gradient sigma
-
bool enabled¶
-
struct PeakFilterParameters¶
- #include <PeakFilter.h>
Parameters for the different filter types.
Public Members
-
double d_min = 1.5¶
minimum d (Bragg’s law)
-
double d_max = 50.0¶
maximum d (Bragg’s law)
-
double strength_min = 1.0¶
minimum strength (I/sigma)
-
double strength_max = 1.0e7¶
maximum strength (I/sigma)
-
double intensity_min = 0.0¶
minimum intensity
-
double intensity_max = 1.0e7¶
maximum intensity
-
double sigma_min = 0.0¶
minimum sigma
-
double sigma_max = 1000.0¶
maximum sigma
-
std::string unit_cell = ""¶
unit cell name
-
double unit_cell_tolerance = 0.2¶
indexing tolerance
-
double significance = 0.99¶
signficance
-
double sparse = 100¶
number of peaks in dataset to be kept
-
double first_frame = 0.0¶
start of frame range
-
double last_frame = 9.0¶
end of frame range
-
double peak_end = 3.0¶
scale for peak intensity ellipsoid (sigmas)
-
double bkg_end = 6.0¶
scale for background ellipsoid (sigmas)
-
RejectionFlag rejection_flag = RejectionFlag::NotRejected¶
rejection flag to keep
-
double d_min = 1.5¶
-
class PeakFilter¶
- #include <PeakFilter.h>
Remove peaks that meet specific criteria from a collection.
Peaks excluded by the filter will be labelled “caught” by the filter.
Public Functions
-
PeakFilter()¶
Constructor.
-
void resetFilterFlags()¶
set filter parameters to default
-
void filterComplementary(PeakCollection *peak_collection) const¶
Filter peaks that are complementary to the given peaks.
-
void filterEnabled(PeakCollection *peak_collection) const¶
Filter only enabled peaks.
-
void filterMasked(PeakCollection *peak_collection) const¶
Filter only masked peaks.
-
void filterIndexed(PeakCollection *peak_collection) const¶
Filter only peaks indexed by their unit cell and its indexing tolerance.
-
void filterIndexTolerance(PeakCollection *peak_collection) const¶
Filter only peaks indexed by a given unit cell with a given indexing tolerance
-
void filterUnitCell(PeakCollection *peak_collection) const¶
Keeps only the peaks whose unit cell is the given unit cell.
-
void filterHasUnitCell(PeakCollection *peak_collection) const¶
Remove peaks whitout unit cell.
-
void filterStrength(PeakCollection *peak_collection) const¶
Filter peaks with I/sigma above threshold.
-
void filterDRange(PeakCollection *peak_collection) const¶
Remove peaks which are not in a d-range.
-
void filterSignificance(PeakCollection *peak_collection) const¶
Filter merged peaks which satisfies a chi2 test.
-
void filterOverlapping(PeakCollection *peak_collection) const¶
Remove overlapping peaks.
-
void filterExtinct(PeakCollection *peak_collection) const¶
Remove space-group extincted peaks.
-
void filterSparseDataSet(PeakCollection *peak_collection) const¶
Remove peaks which belongs to datasets containing too few peaks.
-
void filterFrameRange(PeakCollection *peak_collection) const¶
Remove peaks outside the specified frame range.
-
void filterRejectionFlag(PeakCollection *peak_collection) const¶
Remove peaks without the specified rejection flag.
-
void filterIntensity(PeakCollection *peak_collection) const¶
Remove peaks outside the given intensity range.
-
void filterSigma(PeakCollection *peak_collection) const¶
Remove peaks outside the given sigma range.
-
void filterGradient(PeakCollection *peak_collection) const¶
Remove peaks outside given background gradient range.
-
void filterGradientSigma(PeakCollection *peak_collection) const¶
Remove peaks outside given background gradient sigma range.
-
std::vector<Peak3D*> filterEnabled(const std::vector<Peak3D*> &peaks, bool flag) const¶
Filter only enabled on a peak vector.
-
std::vector<Peak3D*> filterIndexed(const std::vector<Peak3D*> &peaks, const UnitCell *cell = nullptr, const InstrumentState *state = nullptr) const¶
Filter only enabled on a peak vector.
-
std::vector<Peak3D*> filterFrameRange(const std::vector<Peak3D*> &peaks, int first_frame, int last_frame) const¶
Remove peaks from outside a frame range for a vector.
-
std::vector<Peak3D*> filterDRange(const std::vector<Peak3D*> &peaks, double d_min, double d_max, const InstrumentState *state = nullptr) const¶
Filter d-range on peak vector (intended for a single frame, hence the state)
-
std::vector<Peak3D*> filterStrength(const std::vector<Peak3D*> &peaks, double str_min, double str_max)¶
Filter strength on peak vector.
-
void filter(PeakCollection *peak_collection) const¶
Run the filtering.
-
void resetFiltering(PeakCollection *peak_collection) const¶
Reset filtering.
-
PeakFilterParameters *parameters()¶
Get a pointer to the filter parameters.
-
PeakFilterFlags *flags()¶
Get pointer to filter flags.
-
PeakFilter()¶
-
struct PredictionParameters : public ohkl::IntegrationParameters¶
- #include <Predictor.h>
Parameters for peak prediction.
Subclassed by ohkl::StrategyParameters
Public Members
-
double d_min = 1.5¶
Minimum detector range (filter)
-
double d_max = 50.0¶
Maximum detector range (filter)
-
double peak_end = 3.0¶
End of peak region for RegionType::VariableEllipsoid (sigmas)
-
double bkg_begin = 3.0¶
Beginning of background region (sigmas)
-
double bkg_end = 6.0¶
End of background region (sigmas)
-
double fixed_peak_end = 8¶
End of peak region for RegionType::FixedEllipsoid (pixels)
-
double fixed_bkg_begin = 1.0¶
Beginning of background region (factor of fixed_peak_end)
-
double fixed_bkg_end = 2.0¶
End of background region (factor of fixed_peak_end)
-
double max_counts = 50000.0¶
Maximum per-pixel count.
-
double max_strength = 1¶
Maximum strength for profile integration.
-
double max_d = 2.5¶
Maximum d (minimum resolution) for profile integration.
-
int max_width = 10¶
Maximum number of images a peak/shape can span.
-
bool fit_center = true¶
Whether to update the peak centre after integration.
-
bool fit_cov = true¶
Whether to update the peak covariance after integration.
-
bool discard_saturated = false¶
Whether to discard peaks with saturated pixels.
-
IntegratorType integrator_type = IntegratorType::PixelSum¶
Type of integrator.
-
bool use_gradient = false¶
Use gradient to discriminate heterogeneous background regions.
-
GradientKernel gradient_type = GradientKernel::Sobel¶
Kernel to use for gradient convolution.
-
bool fft_gradient = false¶
Whether to use FFT or real space gradient computation.
-
RegionType region_type = RegionType::VariableEllipsoid¶
Whether to use fixed or sigma-dependent integration regions.
-
bool skip_masked = true¶
Whether to skip peaks that intersect masks.
-
bool remove_overlaps = false¶
Whether to remove peaks with overlapping peak integration regions.
-
bool use_max_strength = false¶
Skip profile integration of peaks with strength above threshold.
-
bool use_max_d = false¶
Skip profile integration of peaks with resolution below threshold.
-
bool use_max_width = false¶
Skip integration of peaks spanning too many images.
-
double d_min = 1.5¶
-
struct StrategyParameters : public ohkl::PredictionParameters¶
- #include <Predictor.h>
Parameters for strategy tool.
Public Members
-
double d_min = 1.5¶
Minimum detector range (filter)
-
double d_max = 50.0¶
Maximum detector range (filter)
-
double peak_end = 3.0¶
End of peak region for RegionType::VariableEllipsoid (sigmas)
-
double bkg_begin = 3.0¶
Beginning of background region (sigmas)
-
double bkg_end = 6.0¶
End of background region (sigmas)
-
double fixed_peak_end = 8¶
End of peak region for RegionType::FixedEllipsoid (pixels)
-
double fixed_bkg_begin = 1.0¶
Beginning of background region (factor of fixed_peak_end)
-
double fixed_bkg_end = 2.0¶
End of background region (factor of fixed_peak_end)
-
double max_counts = 50000.0¶
Maximum per-pixel count.
-
double max_strength = 1¶
Maximum strength for profile integration.
-
double max_d = 2.5¶
Maximum d (minimum resolution) for profile integration.
-
int max_width = 10¶
Maximum number of images a peak/shape can span.
-
bool fit_center = true¶
Whether to update the peak centre after integration.
-
bool fit_cov = true¶
Whether to update the peak covariance after integration.
-
bool discard_saturated = false¶
Whether to discard peaks with saturated pixels.
-
IntegratorType integrator_type = IntegratorType::PixelSum¶
Type of integrator.
-
bool use_gradient = false¶
Use gradient to discriminate heterogeneous background regions.
-
GradientKernel gradient_type = GradientKernel::Sobel¶
Kernel to use for gradient convolution.
-
bool fft_gradient = false¶
Whether to use FFT or real space gradient computation.
-
RegionType region_type = RegionType::VariableEllipsoid¶
Whether to use fixed or sigma-dependent integration regions.
-
bool skip_masked = true¶
Whether to skip peaks that intersect masks.
-
bool remove_overlaps = false¶
Whether to remove peaks with overlapping peak integration regions.
-
bool use_max_strength = false¶
Skip profile integration of peaks with strength above threshold.
-
bool use_max_d = false¶
Skip profile integration of peaks with resolution below threshold.
-
bool use_max_width = false¶
Skip integration of peaks spanning too many images.
-
double d_min = 1.5¶
-
class Predictor¶
- #include <Predictor.h>
Predict peaks positions in real space given a unit cell.
Given a unit cell, convert all combinations of Miller indices in some range to q-vectors, and transform their positions to real space detector spots.
Public Functions
-
void predictPeaks(const sptrDataSet data, const sptrUnitCell unit_cell)¶
Predict peaks give a unit cell.
-
std::vector<Peak3D*> buildPeaksFromMillerIndices(sptrDataSet data, const std::vector<MillerIndex> &hkls, const sptrUnitCell unit_cell, sptrProgressHandler handler = nullptr)¶
Build a list of peaks from hkls as computed from unit cell.
-
void strategyPredict(sptrDataSet data, const sptrUnitCell unit_cell)¶
Predict peaks in strategy mode.
-
PredictionParameters *parameters()¶
Get a pointer to the prediction parameters.
-
StrategyParameters *strategyParamters()¶
Get a pointer to the strategy paramters.
-
unsigned int numberOfPredictedPeaks()¶
Get the number of predicted peaks.
-
void setHandler(sptrProgressHandler handler)¶
Set handler for GUI.
-
void predictPeaks(const sptrDataSet data, const sptrUnitCell unit_cell)¶
-
class Profile1D¶
- #include <Profile1D.h>
A 1D peak profile used for profile-fitting integration.
A profile is the average of many peaks. This class allows profiles to be summed and normalised for use in integration.
Public Functions
-
Profile1D(const Intensity &mean_background = Intensity(), double sigma_max = 4.0, size_t num = 200)¶
Constructor. sigma_max indicates maximum number of standard deviations.
-
void addPoint(double r2, double M)¶
Add a data point to the bins.
- Parameters:
r2 – : equal to \((x-x0)\cdot A\cdot (x-x0)\) where \(A\) is the inverse covariance matrix of the peak.
M – : the total count (no background correction)
-
const std::vector<double> &counts() const¶
Returns the vector of integrated counts values.
-
const std::vector<int> &npoints() const¶
Returns the number of points in each bin.
-
std::vector<Intensity> profile() const¶
Returns the profile \(I(s) / I(s_{max})\).
-
void reset()¶
Reset the profile to zero.
-
Profile1D(const Intensity &mean_background = Intensity(), double sigma_max = 4.0, size_t num = 200)¶
-
class Profile3D¶
- #include <Profile3D.h>
A 3D peak profile used for profile-fitting integration.
A profile is the average of many peaks. This class allows profiles to be summed and normalised for use in integration.
Public Functions
-
Profile3D(const AABB &aabb, int nx, int ny, int nz)¶
Construct with given bounding box and number of bins
- Parameters:
aabb – the bounding box for the profile
nx – : number of histogram bins in x direction
ny – : number of histogram bins in y direction
nz – : number of histogram bins in z (frame) direction
-
Profile3D(const AABB &aabb, const Eigen::Vector3i &shape)¶
Construct with given bounding box and shape tensor.
-
double at(size_t i, size_t j, size_t k) const¶
Return the value of the bin at the given indices.
-
double &operator()(size_t i, size_t j, size_t k)¶
Return the value of the bin at the given indices.
-
const double &operator()(size_t i, size_t j, size_t k) const¶
Return the value of the bin at the given indices.
-
bool addValue(const Eigen::Vector3d &x, double y)¶
Add a value to the profile.
-
size_t addSubdividedValue(const Eigen::Vector3d &x, double y, size_t subdivide)¶
Add a value to the profile, with subdivision.
-
const Eigen::Vector3d &dx() const¶
Compute the bin size.
-
size_t count() const¶
Returns the number of points in the profile.
-
double predict(const Eigen::Vector3d &x) const¶
Predict the count value for a given histogram bin.
-
bool normalize()¶
Normalize the profile so that the sum is equal to one. Returns false if it cannot be normalized.
-
void addProfile(const Profile3D &other, double weight)¶
Add the given profile to this one, with the specified weight.
-
Ellipsoid ellipsoid() const¶
Compute the ellipsoid corresponding to the center of mass and inertia tensor.
-
const AABB &aabb() const¶
Returns the bounding box of the profile.
-
const Eigen::Vector3i &shape() const¶
Returns the shape of the underlying data.
-
inline int nProfiles() const¶
Return the number of profiles used to construct this one (via meanProfile)
-
Profile3D(const AABB &aabb, int nx, int ny, int nz)¶
-
struct ShapeModelParameters : public ohkl::IntegrationParameters¶
- #include <ShapeModel.h>
Parameters for building the shape collection.
Public Functions
-
void log(const Level &level) const¶
Search radius for neighbouring profiles (frames)
Public Members
-
double d_min = 1.5¶
Minimum d (filter)
-
double d_max = 50.0¶
Maximum d (filter)
-
double strength_min = 1.0¶
Minimum peak strength I/sigma (filter)
-
double strength_max = 1.0e7¶
Maximum peak strength I/sigma (filter)
-
bool kabsch_coords = true¶
Are we using Kabsch or detector coordinates?
-
int nbins_x = 20¶
Number of x histogram bins for peak.
-
int nbins_y = 20¶
Number of y histogram bins for peak.
-
int nbins_z = 10¶
Number of z histogram bins for peak.
-
int n_subdiv = 1¶
Number subdivisions along each axis per pixel.
-
int min_n_neighbors = 10¶
Minimum number of neighbours required for shape collection.
-
double sigma_m = 0.1¶
Variance due to crystal mosaicity.
-
double sigma_d = 0.1¶
Variance due to beam divergence.
-
double neighbour_range_frames = 10.0¶
Search radius for neighbouring profiles (pixels)
-
double peak_end = 3.0¶
End of peak region for RegionType::VariableEllipsoid (sigmas)
-
double bkg_begin = 3.0¶
Beginning of background region (sigmas)
-
double bkg_end = 6.0¶
End of background region (sigmas)
-
double fixed_peak_end = 8¶
End of peak region for RegionType::FixedEllipsoid (pixels)
-
double fixed_bkg_begin = 1.0¶
Beginning of background region (factor of fixed_peak_end)
-
double fixed_bkg_end = 2.0¶
End of background region (factor of fixed_peak_end)
-
double max_counts = 50000.0¶
Maximum per-pixel count.
-
double max_strength = 1¶
Maximum strength for profile integration.
-
double max_d = 2.5¶
Maximum d (minimum resolution) for profile integration.
-
int max_width = 10¶
Maximum number of images a peak/shape can span.
-
bool fit_center = true¶
Whether to update the peak centre after integration.
-
bool fit_cov = true¶
Whether to update the peak covariance after integration.
-
bool discard_saturated = false¶
Whether to discard peaks with saturated pixels.
-
IntegratorType integrator_type = IntegratorType::PixelSum¶
Type of integrator.
-
bool use_gradient = false¶
Use gradient to discriminate heterogeneous background regions.
-
GradientKernel gradient_type = GradientKernel::Sobel¶
Kernel to use for gradient convolution.
-
bool fft_gradient = false¶
Whether to use FFT or real space gradient computation.
-
RegionType region_type = RegionType::VariableEllipsoid¶
Whether to use fixed or sigma-dependent integration regions.
-
bool skip_masked = true¶
Whether to skip peaks that intersect masks.
-
bool remove_overlaps = false¶
Whether to remove peaks with overlapping peak integration regions.
-
bool use_max_strength = false¶
Skip profile integration of peaks with strength above threshold.
-
bool use_max_d = false¶
Skip profile integration of peaks with resolution below threshold.
-
bool use_max_width = false¶
Skip integration of peaks spanning too many images.
-
void log(const Level &level) const¶
-
class ShapeModel¶
- #include <ShapeModel.h>
Store a collection of peak shapes, for peak prediction and integration.
The collection stores a list of reference peaks. For each reference peak, the collection stores the covariance matrix of the intensity distribution, as well as 3d- and 1d- integrated profiles. The covariance matrices are used to roughly predict the shapes of weak peaks, and the integrated profiles are used in the profile-fitting integration methods.
Public Functions
-
ShapeModel(const sptrDataSet data)¶
Construct an empty collection.
- Parameters:
detector_coords – if true, store profiles in detector coordinates; otherwise store in Kabsch coordinates
-
inline unsigned int id() const¶
Get the integer id.
-
void setId(unsigned int id)¶
Set the integer id.
-
inline void setName(const std::string &name)¶
Set the name.
-
inline std::string name()¶
Get the name.
-
bool addPeak(Peak3D *peak, Profile3D &&profile, Profile1D &&integrated_profile)¶
Add a reference peak to the collection.
-
void updateFit(int num_iterations = 1000)¶
Update the fitted covariances.
-
void setParameters(const ShapeModelParameters ¶ms)¶
Set the shape collection parameters.
-
void setPredictedShapes(PeakCollection *peaks, bool thread_parallel = true)¶
Set shapes of a predicted peak collection.
-
Eigen::Matrix3d predictCovariance(Peak3D *peak) const¶
Predict the (detector space) covariance of a given peak.
-
double meanPearson() const¶
Returns mean Pearson coefficient to measure quality of fit.
-
Profile3D meanProfile(const DetectorEvent &ev) const¶
Returns the average or nearest peak profile near the given detector event.
-
std::vector<Intensity> meanProfile1D(const DetectorEvent &ev) const¶
Returns the average or nearest peak profile near the given detector event.
-
Eigen::Matrix3d meanCovariance(Peak3D *reference_peak) const¶
Returns the average or nearest peak covariance near the given detector event.
-
std::vector<Peak3D*> findNeighbors(const DetectorEvent &ev) const¶
Find neighbors of a given peak, or the nearest peak if there are none within the cutoff.
-
std::array<double, 6> choleskyD() const¶
Returns the background end used for the collection.
-
std::array<double, 6> choleskyM() const¶
Returns the background end used for the collection.
-
std::array<double, 6> choleskyS() const¶
Returns the background end used for the collection.
-
std::map<Peak3D*, std::pair<Profile3D, Profile1D>> profiles() const¶
Returns the background end used for the collection.
-
inline int numberOfPeaks() const¶
Return number of peaks in collection.
-
ShapeModelParameters *parameters()¶
Shape collection parameters.
-
ShapeModel(const sptrDataSet data)¶
-
using RankedSolution = std::pair<sptrUnitCell, double>¶