Mantid
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
Mantid::MDAlgorithms::Integrate3DEvents Class Reference

This is a low-level class to construct a map with lists of events near each peak Q-vector, shifted to be centered at (0,0,0). More...

#include <Integrate3DEvents.h>

Public Member Functions

void addEvents (std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &event_qs, bool hkl_integ)
 Add event Q's to lists of events near peaks. More...
 
std::shared_ptr< const Mantid::Geometry::PeakShapeellipseIntegrateEvents (const std::vector< Kernel::V3D > &E1Vec, Mantid::Kernel::V3D const &peak_q, bool specify_size, double peak_radius, double back_inner_radius, double back_outer_radius, std::vector< double > &axes_radii, double &inti, double &sigi)
 Find the net integrated intensity of a peak, using ellipsoidal volumes. More...
 
std::shared_ptr< const Mantid::Geometry::PeakShapeellipseIntegrateModEvents (const std::vector< Kernel::V3D > &E1Vec, Mantid::Kernel::V3D const &peak_q, Mantid::Kernel::V3D const &hkl, Mantid::Kernel::V3D const &mnp, bool specify_size, double peak_radius, double back_inner_radius, double back_outer_radius, std::vector< double > &axes_radii, double &inti, double &sigi)
 Find the net integrated intensity of a modulated peak, using ellipsoidal volumes. More...
 
double estimateSignalToNoiseRatio (const IntegrationParameters &params, const Mantid::Kernel::V3D &center, bool forceSpherical=false, double sphericityTol=0.02)
 
 Integrate3DEvents (const std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > &peak_q_list, Kernel::DblMatrix UBinv, double radius, const bool useOnePercentBackgroundCorrection=true)
 Construct object to store events around peaks and integrate peaks. More...
 
 Integrate3DEvents (const std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > &peak_q_list, std::vector< Mantid::Kernel::V3D > const &hkl_list, std::vector< Mantid::Kernel::V3D > const &mnp_list, Kernel::DblMatrix UBinv, Kernel::DblMatrix ModHKL, double radius_m, double radius_s, int MaxO, const bool CrossT, const bool useOnePercentBackgroundCorrection=true)
 With modulation vectors, construct an object to store events that correspond to a peak and are within the specified radius of the specified peak centers, and to integrate the peaks. More...
 
std::pair< std::shared_ptr< const Mantid::Geometry::PeakShape >, std::tuple< double, double, double > > integrateStrongPeak (const IntegrationParameters &params, const Kernel::V3D &peak_q, double &inti, double &sigi)
 Find the net integrated intensity of a peak, using ellipsoidal volumes. More...
 
std::shared_ptr< const Geometry::PeakShapeintegrateWeakPeak (const IntegrationParameters &params, Mantid::DataObjects::PeakShapeEllipsoid_const_sptr shape, const std::tuple< double, double, double > &libPeak, const Mantid::Kernel::V3D &peak_q, double &inti, double &sigi)
 

Private Member Functions

void addEvent (std::pair< std::pair< double, double >, Mantid::Kernel::V3D > event_Q, bool hkl_integ)
 Add an event to the vector of events for the closest h,k,l. More...
 
void addModEvent (std::pair< std::pair< double, double >, Mantid::Kernel::V3D > event_Q, bool hkl_integ)
 Add an event to the appropriate vector of events for the closest h,k,l, if it is within the required radius of the corresponding peak in the PeakQMap. More...
 
std::tuple< double, double, double > calculateRadiusFactors (const IntegrationParameters &params, double max_sigma) const
 Calculate the radius to use for each axis of the ellipsoid from the parameters provided. More...
 
bool correctForDetectorEdges (std::tuple< double, double, double > &radii, const std::vector< Mantid::Kernel::V3D > &E1Vecs, const Mantid::Kernel::V3D &peak_q, const std::vector< double > &axesRadii, const std::vector< double > &bkgInnerRadii, const std::vector< double > &bkgOuterRadii)
 
double detectorQ (const std::vector< Kernel::V3D > &E1Vec, const Mantid::Kernel::V3D &QLabFrame, const std::vector< double > &r)
 Compute if a particular Q falls on the edge of a detector. More...
 
std::shared_ptr< const Mantid::DataObjects::PeakShapeEllipsoidellipseIntegrateEvents (const std::vector< Kernel::V3D > &E1Vec, Kernel::V3D const &peak_q, std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &ev_list, std::vector< Mantid::Kernel::V3D > const &directions, std::vector< double > const &sigmas, bool specify_size, double peak_radius, double back_inner_radius, double back_outer_radius, std::vector< double > &axes_radii, double &inti, double &sigi)
 Find the net integrated intensity of a list of Q's using ellipsoids. More...
 
const std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > * getEvents (const Mantid::Kernel::V3D &peak_q)
 Get a list of events for a given Q. More...
 
int64_t getHklKey (Mantid::Kernel::V3D const &q_vector)
 Form a map key for the specified q_vector. More...
 
int64_t getHklKey2 (Mantid::Kernel::V3D const &hkl)
 Form a map key for the specified q_vector. More...
 
int64_t getHklMnpKey (Mantid::Kernel::V3D const &q_vector)
 Form a map key for the specified q_vector of satellite peaks. More...
 
int64_t getHklMnpKey2 (Mantid::Kernel::V3D const &hkl)
 Form a map key for the specified q_vector. More...
 

Static Private Member Functions

static void getEigenVectors (Kernel::DblMatrix const &cov_matrix, std::vector< Mantid::Kernel::V3D > &eigen_vectors, std::vector< double > &eigen_values)
 Calculate the eigen vectors of a 3x3 real symmetric matrix. More...
 
static int64_t getHklKey (int h, int k, int l)
 Form a map key as 10^12*h + 10^6*k + l from the integers h, k, l. More...
 
static int64_t getHklMnpKey (int h, int k, int l, int m, int n, int p)
 Form a map key as 10^12*h + 10^6*k + l from the integers h,k,l. More...
 
static void makeCovarianceMatrix (std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &events, Kernel::DblMatrix &matrix, double radius)
 Calculate the 3x3 covariance matrix of a list of Q-vectors at 0,0,0. More...
 
static std::pair< double, double > numInEllipsoid (std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &events, std::vector< Mantid::Kernel::V3D > const &directions, std::vector< double > const &sizes)
 Calculate the number of events in an ellipsoid centered at 0,0,0. More...
 
static std::pair< double, double > numInEllipsoidBkg (std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &events, std::vector< Mantid::Kernel::V3D > const &directions, std::vector< double > const &sizes, std::vector< double > const &sizesIn, const bool useOnePercentBackgroundCorrection)
 Calculate the number of events in an ellipsoid centered at 0,0,0. More...
 

Private Attributes

const bool crossterm
 
EventListMap m_event_lists
 
Kernel::DblMatrix m_ModHKL
 
PeakQMap m_peak_qs
 
double m_radius
 
Kernel::DblMatrix m_UBinv
 
const bool m_useOnePercentBackgroundCorrection
 
int maxOrder
 
double s_radius
 

Detailed Description

This is a low-level class to construct a map with lists of events near each peak Q-vector, shifted to be centered at (0,0,0).

A method is also provided to find the principal axes of such a list of events, and to find the net integrated counts, using ellipsoids with axis lengths determined from the standard deviations in the directions of the principal axes.

Author
Dennis Mikkelson
Date
2012-12-19

Definition at line 58 of file Integrate3DEvents.h.

Constructor & Destructor Documentation

◆ Integrate3DEvents() [1/2]

Mantid::MDAlgorithms::Integrate3DEvents::Integrate3DEvents ( const std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > &  peak_q_list,
Kernel::DblMatrix  UBinv,
double  radius,
const bool  useOnePercentBackgroundCorrection = true 
)

Construct object to store events around peaks and integrate peaks.

Construct an object to store events that correspond to a peak and are within the specified radius of the specified peak centers, and to integrate the peaks.

Parameters
peak_q_listList of Q-vectors for peak centers.
UBinvThe matrix that maps Q-vectors to h,k,l
radiusThe maximum distance from a peak's Q-vector, for an event to be stored in the list associated with that peak.
useOnePercentBackgroundCorrectionflag if one percent background correction should be used.

Definition at line 48 of file Integrate3DEvents.cpp.

References getHklKey(), and m_peak_qs.

◆ Integrate3DEvents() [2/2]

Mantid::MDAlgorithms::Integrate3DEvents::Integrate3DEvents ( const std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > &  peak_q_list,
std::vector< Mantid::Kernel::V3D > const &  hkl_list,
std::vector< Mantid::Kernel::V3D > const &  mnp_list,
Kernel::DblMatrix  UBinv,
Kernel::DblMatrix  ModHKL,
double  radius_m,
double  radius_s,
int  MaxO,
const bool  CrossT,
const bool  useOnePercentBackgroundCorrection = true 
)

With modulation vectors, construct an object to store events that correspond to a peak and are within the specified radius of the specified peak centers, and to integrate the peaks.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Parameters
peak_q_listList of Q-vectors for peak centers.
hkl_listThe list of h,k,l
mnp_listThe list of satellite m,n,p
UBinvThe matrix that maps Q-vectors to h,k,l
ModHKLThe Modulation vectors
radius_mThe maximum distance from a peak's Q-vector, for an event to be stored in the list associated with that peak.
radius_sThe maximum distance from a peak's Q-vector, for an event to be stored in the list associated with that satellite peak.
MaxOThe maximum order of satellite peaks.
CrossTSwitch for cross terms of satellites.
useOnePercentBackgroundCorrectionflag if one percent background correction should be used.

Definition at line 82 of file Integrate3DEvents.cpp.

References getHklMnpKey(), and m_peak_qs.

Member Function Documentation

◆ addEvent()

void Mantid::MDAlgorithms::Integrate3DEvents::addEvent ( std::pair< std::pair< double, double >, Mantid::Kernel::V3D event_Q,
bool  hkl_integ 
)
private

Add an event to the vector of events for the closest h,k,l.

Add an event to the appropriate vector of events for the closest h,k,l, if it is within the required radius of the corresponding peak in the PeakQMap.

NOTE: The event passed in may be modified by this method. In particular, if it corresponds to one of the specified peak_qs, the corresponding peak q will be subtracted from the event and the event will be added to that peak's vector in the event_lists map.

Parameters
event_QThe Q-vector for the event that may be added to the event_lists map, if it is close enough to some peak
hkl_integ

Definition at line 947 of file Integrate3DEvents.cpp.

References getHklKey(), getHklKey2(), m_event_lists, m_peak_qs, m_radius, and m_UBinv.

Referenced by addEvents().

◆ addEvents()

void Mantid::MDAlgorithms::Integrate3DEvents::addEvents ( std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &  event_qs,
bool  hkl_integ 
)

Add event Q's to lists of events near peaks.

Add the specified event Q's to lists of events near peaks.

An event is added to at most one list. First the nearest h,k,l for that event Q vector is calculated. If a peak with that h,k,l was specified when this object was constructed and if the distance from the specified event Q to that peak is less than the radius that was specified at construction time, then the event Q vector is added to the list of event Q vectors for that peak. NOTE: The Q-vectors passed in to this method will be shifted by the center Q for it's associated peak, so that the list of Q-vectors for a peak are centered around 0,0,0 and represent offsets in Q from the peak center.

Parameters
event_qsList of event Q vectors to add to lists of Q's associated with peaks.
hkl_integ

Definition at line 116 of file Integrate3DEvents.cpp.

References addEvent(), addModEvent(), and maxOrder.

Referenced by Mantid::MDAlgorithms::IntegrateEllipsoidsTwoStep::qListFromEventWS(), Mantid::MDAlgorithms::IntegrateEllipsoidsV1::qListFromEventWS(), Mantid::MDAlgorithms::IntegrateEllipsoidsTwoStep::qListFromHistoWS(), and Mantid::MDAlgorithms::IntegrateEllipsoidsV1::qListFromHistoWS().

◆ addModEvent()

void Mantid::MDAlgorithms::Integrate3DEvents::addModEvent ( std::pair< std::pair< double, double >, Mantid::Kernel::V3D event_Q,
bool  hkl_integ 
)
private

Add an event to the appropriate vector of events for the closest h,k,l, if it is within the required radius of the corresponding peak in the PeakQMap.

NOTE: The event passed in may be modified by this method. In particular, if it corresponds to one of the specified peak_qs, the corresponding peak q will be subtracted from the event and the event will be added to that peak's vector in the event_lists map.

Parameters
event_QThe Q-vector for the event that may be added to the event_lists map, if it is close enough to some peak
hkl_integ

Definition at line 985 of file Integrate3DEvents.cpp.

References getHklMnpKey(), getHklMnpKey2(), m_event_lists, m_peak_qs, m_radius, m_UBinv, and s_radius.

Referenced by addEvents().

◆ calculateRadiusFactors()

std::tuple< double, double, double > Mantid::MDAlgorithms::Integrate3DEvents::calculateRadiusFactors ( const IntegrationParameters params,
double  max_sigma 
) const
private

Calculate the radius to use for each axis of the ellipsoid from the parameters provided.

Parameters
params:: the integration parameters
max_sigma:: the largest sigma of all axes
Returns
tuple of values representing the radius for each axis.

Definition at line 1167 of file Integrate3DEvents.cpp.

References Mantid::MDAlgorithms::IntegrationParameters::backgroundInnerRadius, Mantid::MDAlgorithms::IntegrationParameters::backgroundOuterRadius, Mantid::MDAlgorithms::IntegrationParameters::peakRadius, Mantid::MDAlgorithms::IntegrationParameters::regionRadius, and Mantid::MDAlgorithms::IntegrationParameters::specifySize.

Referenced by estimateSignalToNoiseRatio(), integrateStrongPeak(), and integrateWeakPeak().

◆ correctForDetectorEdges()

bool Mantid::MDAlgorithms::Integrate3DEvents::correctForDetectorEdges ( std::tuple< double, double, double > &  radii,
const std::vector< Mantid::Kernel::V3D > &  E1Vecs,
const Mantid::Kernel::V3D peak_q,
const std::vector< double > &  axesRadii,
const std::vector< double > &  bkgInnerRadii,
const std::vector< double > &  bkgOuterRadii 
)
private

Definition at line 345 of file Integrate3DEvents.cpp.

References detectorQ().

Referenced by integrateStrongPeak(), and integrateWeakPeak().

◆ detectorQ()

double Mantid::MDAlgorithms::Integrate3DEvents::detectorQ ( const std::vector< Kernel::V3D > &  E1Vec,
const Mantid::Kernel::V3D QLabFrame,
const std::vector< double > &  r 
)
private

Compute if a particular Q falls on the edge of a detector.

Calculate if this Q is on a detector The distance from C to OE is given by dv=C-E*(C.scalar_prod(E)) If dv.norm<integration_radius, one of the detector trajectories on the edge is too close to the peak This method is applied to all masked pixels.

If there are masked pixels trajectories inside an integration volume, the peak must be rejected.

Parameters
E1VecVector of values for calculating edge of detectors
QLabFrameThe Peak center.
rPeak radius.

Definition at line 1148 of file Integrate3DEvents.cpp.

References Mantid::API::E1(), Mantid::Kernel::V3D::norm(), and Mantid::Kernel::V3D::scalar_prod().

Referenced by correctForDetectorEdges(), and ellipseIntegrateEvents().

◆ ellipseIntegrateEvents() [1/2]

PeakShapeEllipsoid_const_sptr Mantid::MDAlgorithms::Integrate3DEvents::ellipseIntegrateEvents ( const std::vector< Kernel::V3D > &  E1Vec,
Kernel::V3D const &  peak_q,
std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &  ev_list,
std::vector< Mantid::Kernel::V3D > const &  directions,
std::vector< double > const &  sigmas,
bool  specify_size,
double  peak_radius,
double  back_inner_radius,
double  back_outer_radius,
std::vector< double > &  axes_radii,
double &  inti,
double &  sigi 
)
private

Find the net integrated intensity of a list of Q's using ellipsoids.

Integrate a list of events, centered about (0,0,0) given the principal axes for the events and the standard deviations in the the directions of the principal axes.

Parameters
E1VecVector of values for calculating edge of detectors
peak_qThe Q-vector for the peak center.
ev_listList of events centered at (0,0,0) for a particular peak.
directionsThe three principal axes of the list of events
sigmasThe standard deviations of the events in the directions of the three principal axes.
specify_sizeIf true the integration will be done using the ellipsoids with major axes determined by the peak, back_inner and back_outer radii parameters. If false, the integration will be done using a peak region with major axis chosen so that it covers +- three standard deviations of the data in each direction. In this case, the background ellipsoidal shell is chosen to have the same VOLUME as the peak ellipsoid, and to use the peak ellipsoid for the inner radius.
peak_radiusSize of half the major axis of the ellipsoidal peak region.
back_inner_radiusSize of half the major axis of the INNER ellipsoidal boundary of the background region
back_outer_radiusSize of half the major axis of the OUTER ellipsoidal boundary of the background region
axes_radiiThe radii used for integration in the directions of the three principal axes.
intiReturns the net integrated intensity
sigiReturns an estimate of the standard deviation of the net integrated intensity

Definition at line 1051 of file Integrate3DEvents.cpp.

References detectorQ(), m_radius, m_useOnePercentBackgroundCorrection, numInEllipsoid(), numInEllipsoidBkg(), and Mantid::Kernel::QLab.

◆ ellipseIntegrateEvents() [2/2]

Mantid::Geometry::PeakShape_const_sptr Mantid::MDAlgorithms::Integrate3DEvents::ellipseIntegrateEvents ( const std::vector< Kernel::V3D > &  E1Vec,
Mantid::Kernel::V3D const &  peak_q,
bool  specify_size,
double  peak_radius,
double  back_inner_radius,
double  back_outer_radius,
std::vector< double > &  axes_radii,
double &  inti,
double &  sigi 
)

Find the net integrated intensity of a peak, using ellipsoidal volumes.

Integrate the events around the specified peak Q-vector.

The principal axes of the events near this Q-vector and the standard deviations in the directions of these principal axes determine ellipsoidal regions for integrating the peak and estimating the background. Alternatively, if peak and background radii are specified, then those will be used for half the major axis length of the ellipsoids, and the other axes of the ellipsoids will be set proportionally, based on the standard deviations.

Parameters
E1VecVector of values for calculating edge of detectors
peak_qThe Q-vector for the peak center.
specify_sizeIf true the integration will be done using the ellipsoids with major axes determined by the peak, back_inner and back_outer radii parameters. If false, the integration will be done using a peak region with major axis chosen so that it covers +- three standard deviations of the data in each direction. In this case, the background ellipsoidal shell is chosen to have the same VOLUME as the peak ellipsoid, and to use the peak ellipsoid for the inner radius.
peak_radiusSize of half the major axis of the ellipsoidal peak region.
back_inner_radiusSize of half the major axis of the INNER ellipsoidal boundary of the background region
back_outer_radiusSize of half the major axis of the OUTER ellipsoidal boundary of the background region
axes_radiiThe radii used for integration in the directions of the three principal axes.
intiReturns the net integrated intensity
sigiReturns an estimate of the standard deviation of the net integrated intensity

Definition at line 413 of file Integrate3DEvents.cpp.

References ellipseIntegrateEvents(), getEigenVectors(), getHklKey(), m_event_lists, m_radius, makeCovarianceMatrix(), and sigma.

Referenced by ellipseIntegrateEvents(), and ellipseIntegrateModEvents().

◆ ellipseIntegrateModEvents()

Mantid::Geometry::PeakShape_const_sptr Mantid::MDAlgorithms::Integrate3DEvents::ellipseIntegrateModEvents ( const std::vector< Kernel::V3D > &  E1Vec,
Mantid::Kernel::V3D const &  peak_q,
Mantid::Kernel::V3D const &  hkl,
Mantid::Kernel::V3D const &  mnp,
bool  specify_size,
double  peak_radius,
double  back_inner_radius,
double  back_outer_radius,
std::vector< double > &  axes_radii,
double &  inti,
double &  sigi 
)

Find the net integrated intensity of a modulated peak, using ellipsoidal volumes.

Definition at line 462 of file Integrate3DEvents.cpp.

References ellipseIntegrateEvents(), getEigenVectors(), getHklMnpKey(), m_event_lists, m_radius, makeCovarianceMatrix(), s_radius, and sigma.

Referenced by Mantid::MDAlgorithms::IntegrateEllipsoidsV1::exec().

◆ estimateSignalToNoiseRatio()

double Mantid::MDAlgorithms::Integrate3DEvents::estimateSignalToNoiseRatio ( const IntegrationParameters params,
const Mantid::Kernel::V3D center,
bool  forceSpherical = false,
double  sphericityTol = 0.02 
)

◆ getEigenVectors()

void Mantid::MDAlgorithms::Integrate3DEvents::getEigenVectors ( Kernel::DblMatrix const &  cov_matrix,
std::vector< Mantid::Kernel::V3D > &  eigen_vectors,
std::vector< double > &  eigen_values 
)
staticprivate

Calculate the eigen vectors of a 3x3 real symmetric matrix.

Calculate the eigen vectors of a 3x3 real symmetric matrix using the GSL.

Parameters
cov_matrix3x3 real symmetric matrix.
eigen_vectorsThe eigen vectors for the matrix are returned in this list.
eigen_values3 eigenvalues of matrix

Definition at line 654 of file Integrate3DEvents.cpp.

Referenced by ellipseIntegrateEvents(), ellipseIntegrateModEvents(), estimateSignalToNoiseRatio(), and integrateStrongPeak().

◆ getEvents()

const std::vector< std::pair< std::pair< double, double >, V3D > > * Mantid::MDAlgorithms::Integrate3DEvents::getEvents ( const Mantid::Kernel::V3D peak_q)
private

Get a list of events for a given Q.

Definition at line 326 of file Integrate3DEvents.cpp.

References getHklKey(), getHklMnpKey(), m_event_lists, and maxOrder.

Referenced by estimateSignalToNoiseRatio(), integrateStrongPeak(), and integrateWeakPeak().

◆ getHklKey() [1/2]

int64_t Mantid::MDAlgorithms::Integrate3DEvents::getHklKey ( int  h,
int  k,
int  l 
)
staticprivate

Form a map key as 10^12*h + 10^6*k + l from the integers h, k, l.

Form a map key as 10^12*h + 10^6*k + l from the integers h,k,l.

Parameters
hThe first Miller index
kThe second Miller index
lThe third Miller index

Definition at line 692 of file Integrate3DEvents.cpp.

Referenced by addEvent(), ellipseIntegrateEvents(), getEvents(), getHklKey(), getHklKey2(), and Integrate3DEvents().

◆ getHklKey() [2/2]

int64_t Mantid::MDAlgorithms::Integrate3DEvents::getHklKey ( Mantid::Kernel::V3D const &  q_vector)
private

Form a map key for the specified q_vector.

The q_vector is mapped to h,k,l by UBinv and the map key is then formed from those rounded h,k,l values.

Parameters
q_vectorThe q_vector to be mapped to h,k,l

Definition at line 830 of file Integrate3DEvents.cpp.

References getHklKey(), and m_UBinv.

◆ getHklKey2()

int64_t Mantid::MDAlgorithms::Integrate3DEvents::getHklKey2 ( Mantid::Kernel::V3D const &  hkl)
private

Form a map key for the specified q_vector.

The q_vector is mapped to h,k,l by UBinv and the map key is then formed from those rounded h,k,l values.

Parameters
hklThe q_vector to be mapped to h,k,l

Definition at line 726 of file Integrate3DEvents.cpp.

References getHklKey().

Referenced by addEvent().

◆ getHklMnpKey() [1/2]

int64_t Mantid::MDAlgorithms::Integrate3DEvents::getHklMnpKey ( int  h,
int  k,
int  l,
int  m,
int  n,
int  p 
)
staticprivate

Form a map key as 10^12*h + 10^6*k + l from the integers h,k,l.

Parameters
hThe first Miller index
kThe second Miller index
lThe third Miller index
mThe first modulation index
nThe second modulation index
pThe third modulation index

Definition at line 711 of file Integrate3DEvents.cpp.

References Mantid::Geometry::m, and n.

Referenced by addModEvent(), ellipseIntegrateModEvents(), getEvents(), getHklMnpKey(), getHklMnpKey2(), and Integrate3DEvents().

◆ getHklMnpKey() [2/2]

int64_t Mantid::MDAlgorithms::Integrate3DEvents::getHklMnpKey ( Mantid::Kernel::V3D const &  q_vector)
private

Form a map key for the specified q_vector of satellite peaks.

The q_vector is mapped to h,k,l by UBinv and the map key is then formed from those rounded h+-0.5,k+-0.5,l+-0.5 or h,k,l depending on the offset of satellite peak from main peak.

Parameters
q_vectorThe q_vector to be mapped to h,k,l

Definition at line 846 of file Integrate3DEvents.cpp.

References crossterm, getHklMnpKey(), Mantid::Geometry::m, m_ModHKL, m_radius, m_UBinv, maxOrder, n, s_radius, and Mantid::Geometry::IndexingUtils::ValidIndex().

◆ getHklMnpKey2()

int64_t Mantid::MDAlgorithms::Integrate3DEvents::getHklMnpKey2 ( Mantid::Kernel::V3D const &  hkl)
private

Form a map key for the specified q_vector.

The q_vector is mapped to h,k,l by UBinv and the map key is then formed from those rounded h,k,l values.

Parameters
hklThe q_vector to be mapped to h,k,l

Definition at line 739 of file Integrate3DEvents.cpp.

References crossterm, getHklMnpKey(), Mantid::Geometry::m, m_ModHKL, m_radius, maxOrder, n, s_radius, and Mantid::Geometry::IndexingUtils::ValidIndex().

Referenced by addModEvent().

◆ integrateStrongPeak()

std::pair< std::shared_ptr< const Geometry::PeakShape >, std::tuple< double, double, double > > Mantid::MDAlgorithms::Integrate3DEvents::integrateStrongPeak ( const IntegrationParameters params,
const Kernel::V3D peak_q,
double &  inti,
double &  sigi 
)

◆ integrateWeakPeak()

std::shared_ptr< const Geometry::PeakShape > Mantid::MDAlgorithms::Integrate3DEvents::integrateWeakPeak ( const IntegrationParameters params,
Mantid::DataObjects::PeakShapeEllipsoid_const_sptr  shape,
const std::tuple< double, double, double > &  libPeak,
const Mantid::Kernel::V3D peak_q,
double &  inti,
double &  sigi 
)

◆ makeCovarianceMatrix()

void Mantid::MDAlgorithms::Integrate3DEvents::makeCovarianceMatrix ( std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &  events,
Kernel::DblMatrix matrix,
double  radius 
)
staticprivate

Calculate the 3x3 covariance matrix of a list of Q-vectors at 0,0,0.

Given a list of events, associated with a particular peak and already SHIFTED to be centered at (0,0,0), calculate the 3x3 covariance matrix for finding the principal axes of that local event data.

Only events within the specified radius of (0,0,0) will be used.

The covariance matrix can be easily constructed. X, Y, Z of each peak position are the variables we wish to determine the covariance. The mean position in each dimension has already been calculated on subtracted, since this corresponds to the centre position of each peak, which we knew aprori. The expected values of each correlation test X,X X,Y X,Z e.t.c form the elements of this 3 by 3 matrix, but since the probabilities are equal, we can remove them from the sums of the expected values, and simply divide by the number of events for each matrix element. Note that the diagonal elements form the variance X,X, Y,Y, Z,Z

Parameters
eventsVector of V3D objects containing the Q vectors for a peak, with mean at (0,0,0).
matrixA 3x3 matrix that will be filled out with the covariance matrix for the list of events.
radiusOnly events within this radius of the peak center (0,0,0) will be used for calculating the covariance matrix.

Definition at line 625 of file Integrate3DEvents.cpp.

References radius.

Referenced by ellipseIntegrateEvents(), ellipseIntegrateModEvents(), estimateSignalToNoiseRatio(), and integrateStrongPeak().

◆ numInEllipsoid()

std::pair< double, double > Mantid::MDAlgorithms::Integrate3DEvents::numInEllipsoid ( std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &  events,
std::vector< Mantid::Kernel::V3D > const &  directions,
std::vector< double > const &  sizes 
)
staticprivate

Calculate the number of events in an ellipsoid centered at 0,0,0.

Calculate the number of events in an ellipsoid centered at 0,0,0 with the three specified axes and the three specified sizes in the direction of those axes.

NOTE: The three axes must be mutually orthogonal unit vectors.

Parameters
eventsList of 3D events centered at 0,0,0
directionsList of 3 orthonormal directions for the axes of the ellipsoid.
sizesList of three values a,b,c giving half the length of the three axes of the ellisoid.
Returns
Then number of events that are in or on the specified ellipsoid.

Definition at line 528 of file Integrate3DEvents.cpp.

References count.

Referenced by ellipseIntegrateEvents(), estimateSignalToNoiseRatio(), integrateStrongPeak(), and integrateWeakPeak().

◆ numInEllipsoidBkg()

std::pair< double, double > Mantid::MDAlgorithms::Integrate3DEvents::numInEllipsoidBkg ( std::vector< std::pair< std::pair< double, double >, Mantid::Kernel::V3D > > const &  events,
std::vector< Mantid::Kernel::V3D > const &  directions,
std::vector< double > const &  sizes,
std::vector< double > const &  sizesIn,
const bool  useOnePercentBackgroundCorrection 
)
staticprivate

Calculate the number of events in an ellipsoid centered at 0,0,0.

Calculate the number of events in an ellipsoid centered at 0,0,0 with the three specified axes and the three specified sizes in the direction of those axes.

NOTE: The three axes must be mutually orthogonal unit vectors.

Parameters
eventsList of 3D events centered at 0,0,0
directionsList of 3 orthonormal directions for the axes of the ellipsoid.
sizesList of three values a,b,c giving half the length of the three axes of the ellisoid.
sizesInList of three values a,b,c giving half the length of the three inner axes of the ellisoid.
useOnePercentBackgroundCorrectionflag if one percent background correction should be used.
Returns
Then number of events that are in or on the specified ellipsoid.

Definition at line 564 of file Integrate3DEvents.cpp.

References count.

Referenced by ellipseIntegrateEvents(), estimateSignalToNoiseRatio(), integrateStrongPeak(), and integrateWeakPeak().

Member Data Documentation

◆ crossterm

const bool Mantid::MDAlgorithms::Integrate3DEvents::crossterm
private

Definition at line 168 of file Integrate3DEvents.h.

Referenced by getHklMnpKey(), and getHklMnpKey2().

◆ m_event_lists

EventListMap Mantid::MDAlgorithms::Integrate3DEvents::m_event_lists
private

◆ m_ModHKL

Kernel::DblMatrix Mantid::MDAlgorithms::Integrate3DEvents::m_ModHKL
private

Definition at line 164 of file Integrate3DEvents.h.

Referenced by getHklMnpKey(), and getHklMnpKey2().

◆ m_peak_qs

PeakQMap Mantid::MDAlgorithms::Integrate3DEvents::m_peak_qs
private

Definition at line 161 of file Integrate3DEvents.h.

Referenced by addEvent(), addModEvent(), and Integrate3DEvents().

◆ m_radius

double Mantid::MDAlgorithms::Integrate3DEvents::m_radius
private

◆ m_UBinv

Kernel::DblMatrix Mantid::MDAlgorithms::Integrate3DEvents::m_UBinv
private

Definition at line 163 of file Integrate3DEvents.h.

Referenced by addEvent(), addModEvent(), getHklKey(), and getHklMnpKey().

◆ m_useOnePercentBackgroundCorrection

const bool Mantid::MDAlgorithms::Integrate3DEvents::m_useOnePercentBackgroundCorrection
private
Initial value:
=
true

Definition at line 169 of file Integrate3DEvents.h.

Referenced by ellipseIntegrateEvents(), estimateSignalToNoiseRatio(), integrateStrongPeak(), and integrateWeakPeak().

◆ maxOrder

int Mantid::MDAlgorithms::Integrate3DEvents::maxOrder
private

Definition at line 167 of file Integrate3DEvents.h.

Referenced by addEvents(), getEvents(), getHklMnpKey(), and getHklMnpKey2().

◆ s_radius

double Mantid::MDAlgorithms::Integrate3DEvents::s_radius
private

The documentation for this class was generated from the following files: