BornAgain
1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
|
Abstract detector interface.
Handles also "region of interest" (ROI). In general, the ROI is the whole detector, and all methods related to ROI work on the whole detector. If a ROI different to the detector size is explicitly set, then ROI-related methods work on this reduced ROI. Therefore, when calling ROI related methods, the distinction between "explicit ROI exists: yes/no" does not have to be made by the caller, but it is handled in here. For access to the whole detector, even if a smaller ROI is explicitly defined, use the non-ROI-related methods like totalSize() or axes(). Any method which is not speaking of "explicit ROI" handles the "implicit ROI", i.e. uses an explicitly set ROI or the whole detector if no explicit ROI exists. To access the explicitly set ROI, use the methods which have the explicit in its name, like sizeOfExplicitRegionOfInterest().
Definition at line 57 of file IDetector.h.
Classes | |
struct | RoiOfAxis |
Keeps RegionOfInterest (ROI) data of one axis. More... | |
Public Types | |
using | const_iterator = const SimulationAreaIterator & |
Public Member Functions | |
IDetector () | |
~IDetector () override | |
std::vector< size_t > | active_indices () const |
Returns vector of unmasked detector indices. More... | |
void | addDetAxis (const IAxis &axis) |
void | addMask (const IShape2D &shape, bool mask_value=true) |
Adds mask of given shape to the stack of detector masks. The mask value 'true' means that the channel will be excluded from the simulation. The mask which is added last has priority. More... | |
const PolFilter & | analyzer () const |
Returns detection properties. More... | |
void | applyDetectorResolution (Datafield *p_intensity_map) const |
Applies the detector resolution to the given intensity maps. More... | |
OwningVector< IAxis > | axesClippedToRegionOfInterest () const |
Returns the axes clipped to the region of interest. If no region of interest is explicitly defined, then the whole detector is taken as "region of interest". More... | |
const IAxis & | axis (size_t index) const |
One axis of the complete detector. Any region of interest is not taken into account. More... | |
size_t | axisBinIndex (size_t index, size_t selected_axis) const |
Calculate axis index for given global index. More... | |
SimulationAreaIterator | beginNonMaskedPoints () const |
Create begin-iterator to iterate over all points which are not masked and lay within the "Region of Interest". More... | |
SimulationAreaIterator | beginRegionOfInterestPoints () const |
Create begin-iterator to iterate over all points which lay within the "Region of Interest". If no region of interest is explicitly defined, then the whole detector is taken as "region of interest". No matter whether masked or not. More... | |
void | checkNodeArgs () const |
Raises exception if a parameter value is invalid. More... | |
virtual std::string | className () const =0 |
Returns the class name, to be hard-coded in each leaf class that inherits from INode. More... | |
IDetector * | clone () const override=0 |
std::unique_ptr< DetectorContext > | createContext () const |
Datafield * | createDetectorIntensity (const std::vector< std::unique_ptr< DiffuseElement >> &elements) const |
Returns new intensity map with resolution applied, and cropped to ROI if applicable. More... | |
std::unique_ptr< Datafield > | createDetectorMap () const |
Returns empty detector map in given axes units. This map is a data array limited to the size of the "Region of interest". More... | |
virtual IPixel * | createPixel (size_t index) const =0 |
Creates an IPixel for the given Datafield object and index. More... | |
virtual Coords | defaultCoords () const =0 |
Return default axes units. More... | |
size_t | detectorIndexToRegionOfInterestIndex (size_t detectorIndex) const |
const DetectorMask * | detectorMask () const |
const IDetectorResolution * | detectorResolution () const |
Returns a pointer to detector resolution object. More... | |
SimulationAreaIterator | endNonMaskedPoints () const |
Create end-iterator to iterate over all points which are not masked and lay within the "Region of Interest". More... | |
SimulationAreaIterator | endRegionOfInterestPoints () const |
Create end-iterator to iterate over all points which lay within the "Region of Interest". If no region of interest is explicitly defined, then the whole detector is taken as "region of interest". No matter whether masked or not. More... | |
bool | hasExplicitRegionOfInterest () const |
True if a region of interest is explicitly set. More... | |
virtual size_t | indexOfSpecular (const Beam &beam) const =0 |
Returns index of pixel that contains the specular wavevector. If no pixel contains this specular wavevector, the number of pixels is returned. This corresponds to an overflow index. More... | |
void | iterateOverNonMaskedPoints (std::function< void(const_iterator)> func) const |
Iterate over all non-masked points within "region of interest". If no region of interest is explicitly defined, then the whole detector is taken as "region of interest". More... | |
void | iterateOverRegionOfInterest (std::function< void(const_iterator)> func) const |
Iterate over all points within "region of interest", no matter whether they are masked or not. If no region of interest is explicitly defined, then the whole detector is taken as "region of interest". More... | |
void | maskAll () |
Put the mask for all detector channels (i.e. exclude whole detector from the analysis) More... | |
std::vector< const INode * > | nodeChildren () const override |
Returns all children. More... | |
std::vector< const INode * > | nodeOffspring () const |
Returns all descendants. More... | |
size_t | numberOfElements () const |
Returns number of simulation elements. More... | |
virtual ICoordSystem * | offspecCoords (IAxis *beamAxis, const Direction &beamDirection) const =0 |
virtual std::vector< ParaMeta > | parDefs () const |
Returns the parameter definitions, to be hard-coded in each leaf class. More... | |
size_t | rank () const |
Returns number of defined axes. More... | |
std::pair< double, double > | regionOfInterestBounds (size_t iAxis) const |
The lower and upper bound of the region of interest. If no region of interest is explicitly defined, then the whole detector is taken as "region of interest". More... | |
size_t | regionOfInterestIndexToDetectorIndex (size_t regionOfInterestIndex) const |
Convert an index of the region of interest to an index of the detector. If no region of interest is set, then the index stays unmodified (since ROI == detector area). More... | |
void | resetRegionOfInterest () |
Resets region of interest making whole detector plane available for the simulation. More... | |
virtual CoordSystem2D * | scatteringCoords (const Beam &beam) const =0 |
void | setAnalyzer (R3 direction, double efficiency, double total_transmission) |
Sets the polarization analyzer characteristics of the detector. More... | |
virtual void | setDetectorNormal (const Direction &) |
Inits detector with the beam settings. More... | |
void | setDetectorParameters (size_t n_x, double x_min, double x_max, size_t n_y, double y_min, double y_max) |
Sets equidistant axes. More... | |
void | setDetectorResolution (const IDetectorResolution &p_detector_resolution) |
Sets the detector resolution. More... | |
void | setRegionOfInterest (double xlow, double ylow, double xup, double yup) |
Sets rectangular region of interest with lower left and upper right corners defined. More... | |
void | setResolutionFunction (const IResolutionFunction2D &resFunc) |
size_t | sizeOfRegionOfInterest () const |
The size of the "Region of Interest". Same as totalSize() if no region of interest has been explicitly set. More... | |
size_t | totalSize () const |
Returns total number of pixels. Any region of interest is not taken into account. More... | |
virtual void | transferToCPP () |
Used for Python overriding of clone (see swig/tweaks.py) More... | |
Protected Member Functions | |
IDetector (const IDetector &other) | |
virtual std::string | axisName (size_t index) const =0 |
Returns the name for the axis with given index. More... | |
virtual std::pair< double, double > | boundsOfExplicitRegionOfInterest (size_t iAxis) const |
Lower and upper bound of one axis of an explicitly set ROI. Return 0/0 if no ROI has been explicitly set. More... | |
void | clear () |
size_t | getGlobalIndex (size_t x, size_t y) const |
Calculate global index from two axis indices. More... | |
virtual size_t | sizeOfExplicitRegionOfInterest () const |
Return 0 if no ROI has been explicitly set. Size means number of data points. More... | |
Protected Attributes | |
std::vector< RoiOfAxis > | m_explicitROI |
an explicitly defined region of interest. Empty if no ROI has been defined. Vector index corresponds to axis index in m_axes More... | |
std::vector< double > | m_P |
Private Attributes | |
OwningVector< IAxis > | m_axes |
std::shared_ptr< DetectorMask > | m_detector_mask |
std::unique_ptr< IDetectorResolution > | m_detector_resolution |
PolFilter | m_polAnalyzer |
using IDetector::const_iterator = const SimulationAreaIterator& |
Definition at line 59 of file IDetector.h.
IDetector::IDetector | ( | ) |
Definition at line 42 of file IDetector.cpp.
|
overridedefault |
|
protected |
Definition at line 44 of file IDetector.cpp.
References m_detector_resolution, and setDetectorResolution().
std::vector< size_t > IDetector::active_indices | ( | ) | const |
Returns vector of unmasked detector indices.
Definition at line 346 of file IDetector.cpp.
References SimulationAreaIterator::detectorIndex(), and iterateOverNonMaskedPoints().
Referenced by DetectorContext::setup_context().
void IDetector::addDetAxis | ( | const IAxis & | axis | ) |
Definition at line 57 of file IDetector.cpp.
References axis(), IAxis::clone(), OwningVector< T >::emplace_back(), m_axes, m_detector_mask, and rank().
Referenced by setDetectorParameters().
void IDetector::addMask | ( | const IShape2D & | shape, |
bool | mask_value = true |
||
) |
Adds mask of given shape to the stack of detector masks. The mask value 'true' means that the channel will be excluded from the simulation. The mask which is added last has priority.
Definition at line 360 of file IDetector.cpp.
References m_detector_mask.
Referenced by ISimulation2D::addMask(), and maskAll().
|
inline |
Returns detection properties.
Definition at line 204 of file IDetector.h.
References m_polAnalyzer.
Referenced by ISimulation2D::force_polarized(), ISimulation2D::generateElements(), and DetectorContext::setup_context().
void IDetector::applyDetectorResolution | ( | Datafield * | p_intensity_map | ) | const |
Applies the detector resolution to the given intensity maps.
Definition at line 169 of file IDetector.cpp.
References ASSERT, Frame::cloned_axes(), detectorMask(), Datafield::frame(), iterateOverNonMaskedPoints(), m_detector_resolution, SimulationAreaIterator::roiIndex(), and Datafield::setVector().
Referenced by createDetectorIntensity(), and OffspecSimulation::transferDetectorImage().
OwningVector< IAxis > IDetector::axesClippedToRegionOfInterest | ( | ) | const |
Returns the axes clipped to the region of interest. If no region of interest is explicitly defined, then the whole detector is taken as "region of interest".
Definition at line 136 of file IDetector.cpp.
References axis(), IAxis::clip(), OwningVector< T >::emplace_back(), m_axes, regionOfInterestBounds(), and OwningVector< T >::size().
Referenced by RectangularDetector::offspecCoords(), SphericalDetector::offspecCoords(), RectangularDetector::scatteringCoords(), and SphericalDetector::scatteringCoords().
const IAxis & IDetector::axis | ( | size_t | index | ) | const |
One axis of the complete detector. Any region of interest is not taken into account.
Definition at line 74 of file IDetector.cpp.
References ASSERT, m_axes, and rank().
Referenced by IDetector::RoiOfAxis::RoiOfAxis(), addDetAxis(), axesClippedToRegionOfInterest(), createDetectorMap(), RectangularDetector::createPixel(), SphericalDetector::createPixel(), getGlobalIndex(), RectangularDetector::height(), RectangularDetector::indexOfSpecular(), SphericalDetector::indexOfSpecular(), OffspecSimulation::intensityMapSize(), OffspecSimulation::pack_result(), setRegionOfInterest(), OffspecSimulation::transferDetectorImage(), RectangularDetector::width(), RectangularDetector::xSize(), and RectangularDetector::ySize().
size_t IDetector::axisBinIndex | ( | size_t | index, |
size_t | selected_axis | ||
) | const |
Calculate axis index for given global index.
Definition at line 80 of file IDetector.cpp.
References ASSERT, m_axes, rank(), and OwningVector< T >::size().
Referenced by RectangularDetector::createPixel(), and SphericalDetector::createPixel().
|
protectedpure virtual |
Returns the name for the axis with given index.
Implemented in SphericalDetector, and RectangularDetector.
Referenced by setDetectorParameters().
SimulationAreaIterator IDetector::beginNonMaskedPoints | ( | ) | const |
Create begin-iterator to iterate over all points which are not masked and lay within the "Region of Interest".
Definition at line 261 of file IDetector.cpp.
References SimulationAreaIterator::createBegin(), and SimulationAreaIterator::notMasked.
Referenced by createDetectorIntensity(), and iterateOverNonMaskedPoints().
SimulationAreaIterator IDetector::beginRegionOfInterestPoints | ( | ) | const |
Create begin-iterator to iterate over all points which lay within the "Region of Interest". If no region of interest is explicitly defined, then the whole detector is taken as "region of interest". No matter whether masked or not.
Definition at line 271 of file IDetector.cpp.
References SimulationAreaIterator::createBegin(), and SimulationAreaIterator::regionOfInterest.
Referenced by iterateOverRegionOfInterest().
|
protectedvirtual |
Lower and upper bound of one axis of an explicitly set ROI. Return 0/0 if no ROI has been explicitly set.
Definition at line 106 of file IDetector.cpp.
References ASSERT, m_explicitROI, and rank().
Referenced by regionOfInterestBounds().
|
inherited |
Raises exception if a parameter value is invalid.
Definition at line 27 of file INode.cpp.
References ASSERT, RealLimits::check(), INode::className(), INF, RealLimits::limited(), RealLimits::limitless(), INode::m_P, ParaMeta::name, RealLimits::nonnegative(), INode::parDefs(), ParaMeta::vMax, and ParaMeta::vMin.
Referenced by BarGauss::BarGauss(), BarLorentz::BarLorentz(), Bipyramid4::Bipyramid4(), Box::Box(), CantellatedCube::CantellatedCube(), Cone::Cone(), ConstantBackground::ConstantBackground(), CosineRippleBox::CosineRippleBox(), CosineRippleGauss::CosineRippleGauss(), CosineRippleLorentz::CosineRippleLorentz(), Cylinder::Cylinder(), DistributionCosine::DistributionCosine(), DistributionGate::DistributionGate(), DistributionGaussian::DistributionGaussian(), DistributionLogNormal::DistributionLogNormal(), DistributionLorentz::DistributionLorentz(), DistributionTrapezoid::DistributionTrapezoid(), Dodecahedron::Dodecahedron(), EllipsoidalCylinder::EllipsoidalCylinder(), FootprintGauss::FootprintGauss(), FootprintSquare::FootprintSquare(), FuzzySphere::FuzzySphere(), GaussSphere::GaussSphere(), HemiEllipsoid::HemiEllipsoid(), HollowSphere::HollowSphere(), HorizontalCylinder::HorizontalCylinder(), Icosahedron::Icosahedron(), LongBoxGauss::LongBoxGauss(), LongBoxLorentz::LongBoxLorentz(), PlatonicOctahedron::PlatonicOctahedron(), PlatonicTetrahedron::PlatonicTetrahedron(), Prism3::Prism3(), Prism6::Prism6(), Profile1DCauchy::Profile1DCauchy(), Profile1DCosine::Profile1DCosine(), Profile1DGate::Profile1DGate(), Profile1DGauss::Profile1DGauss(), Profile1DTriangle::Profile1DTriangle(), Profile1DVoigt::Profile1DVoigt(), Profile2DCauchy::Profile2DCauchy(), Profile2DCone::Profile2DCone(), Profile2DGate::Profile2DGate(), Profile2DGauss::Profile2DGauss(), Profile2DVoigt::Profile2DVoigt(), Pyramid2::Pyramid2(), Pyramid3::Pyramid3(), Pyramid4::Pyramid4(), Pyramid6::Pyramid6(), RotationEuler::RotationEuler(), RotationX::RotationX(), RotationY::RotationY(), RotationZ::RotationZ(), SawtoothRippleBox::SawtoothRippleBox(), SawtoothRippleGauss::SawtoothRippleGauss(), SawtoothRippleLorentz::SawtoothRippleLorentz(), Sphere::Sphere(), Spheroid::Spheroid(), TruncatedCube::TruncatedCube(), TruncatedSphere::TruncatedSphere(), and TruncatedSpheroid::TruncatedSpheroid().
|
pure virtualinherited |
Returns the class name, to be hard-coded in each leaf class that inherits from INode.
Implemented in Particle, SphericalDetector, SpecularSimulation, ScatteringSimulation, OffspecSimulation, DepthProbeSimulation, PoissonBackground, ConstantBackground, GaussSphere, FuzzySphere, RotationEuler, RotationZ, RotationY, RotationX, IdentityRotation, ParticleCoreShell, ParticleComposition, MesoCrystal, Crystal, MultiLayer, Layer, Lattice3D, HexagonalLattice2D, SquareLattice2D, BasicLattice2D, LayerRoughness, LayerInterface, TruncatedSpheroid, TruncatedSphere, TruncatedCube, Spheroid, Sphere, SawtoothRippleLorentz, SawtoothRippleGauss, SawtoothRippleBox, Pyramid6, Pyramid4, Pyramid3, Pyramid2, Prism6, Prism3, PlatonicTetrahedron, PlatonicOctahedron, LongBoxLorentz, LongBoxGauss, Icosahedron, HorizontalCylinder, HollowSphere, HemiEllipsoid, EllipsoidalCylinder, Dodecahedron, Cylinder, CosineRippleLorentz, CosineRippleGauss, CosineRippleBox, Cone, CantellatedCube, Box, Bipyramid4, BarLorentz, BarGauss, Profile2DVoigt, Profile2DCone, Profile2DGate, Profile2DGauss, Profile2DCauchy, Profile1DVoigt, Profile1DCosine, Profile1DTriangle, Profile1DGate, Profile1DGauss, Profile1DCauchy, MisesGaussPeakShape, MisesFisherGaussPeakShape, LorentzFisherPeakShape, GaussFisherPeakShape, IsotropicLorentzPeakShape, IsotropicGaussPeakShape, ParticleLayout, InterferenceTwin, InterferenceRadialParaCrystal, InterferenceNone, InterferenceHardDisk, InterferenceFinite3DLattice, InterferenceFinite2DLattice, Interference3DLattice, Interference2DSuperLattice, Interference2DParaCrystal, Interference2DLattice, Interference1DLattice, DistributionTrapezoid, DistributionCosine, DistributionLogNormal, DistributionGaussian, DistributionLorentz, DistributionGate, ResolutionFunction2DGaussian, ConvolutionDetectorResolution, PolFilter, RectangularDetector, FootprintSquare, FootprintGauss, and Beam.
Referenced by INode::checkNodeArgs(), ExemplarySamples::createBasic2DParaCrystalWithFTDis(), IProfile1D::pythonConstructor(), IProfile2D::pythonConstructor(), IFormFactor::pythonConstructor(), and IFormFactor::shapeName().
|
protected |
Definition at line 69 of file IDetector.cpp.
References OwningVector< T >::clear(), and m_axes.
Referenced by setDetectorParameters().
|
overridepure virtual |
Implements ICloneable.
Implemented in SphericalDetector, and RectangularDetector.
std::unique_ptr< DetectorContext > IDetector::createContext | ( | ) | const |
Definition at line 355 of file IDetector.cpp.
Referenced by ISimulation2D::prepareSimulation().
Datafield * IDetector::createDetectorIntensity | ( | const std::vector< std::unique_ptr< DiffuseElement >> & | elements | ) | const |
Returns new intensity map with resolution applied, and cropped to ROI if applicable.
Definition at line 191 of file IDetector.cpp.
References applyDetectorResolution(), ASSERT, beginNonMaskedPoints(), createDetectorMap(), endNonMaskedPoints(), and m_detector_resolution.
std::unique_ptr< Datafield > IDetector::createDetectorMap | ( | ) | const |
Returns empty detector map in given axes units. This map is a data array limited to the size of the "Region of interest".
Definition at line 205 of file IDetector.cpp.
References ASSERT, axis(), IAxis::clip(), IAxis::clone(), rank(), and regionOfInterestBounds().
Referenced by createDetectorIntensity().
|
pure virtual |
Creates an IPixel for the given Datafield object and index.
Implemented in SphericalDetector, and RectangularDetector.
Referenced by DetectorContext::setup_context().
|
pure virtual |
Return default axes units.
Implemented in SphericalDetector, and RectangularDetector.
size_t IDetector::detectorIndexToRegionOfInterestIndex | ( | size_t | detectorIndex | ) | const |
const DetectorMask * IDetector::detectorMask | ( | ) | const |
Definition at line 372 of file IDetector.cpp.
References m_detector_mask.
Referenced by applyDetectorResolution(), and SimulationAreaIterator::isMasked().
const IDetectorResolution * IDetector::detectorResolution | ( | ) | const |
Returns a pointer to detector resolution object.
Definition at line 186 of file IDetector.cpp.
References m_detector_resolution.
SimulationAreaIterator IDetector::endNonMaskedPoints | ( | ) | const |
Create end-iterator to iterate over all points which are not masked and lay within the "Region of Interest".
Definition at line 266 of file IDetector.cpp.
References SimulationAreaIterator::createEnd(), and SimulationAreaIterator::notMasked.
Referenced by createDetectorIntensity(), and iterateOverNonMaskedPoints().
SimulationAreaIterator IDetector::endRegionOfInterestPoints | ( | ) | const |
Create end-iterator to iterate over all points which lay within the "Region of Interest". If no region of interest is explicitly defined, then the whole detector is taken as "region of interest". No matter whether masked or not.
Definition at line 276 of file IDetector.cpp.
References SimulationAreaIterator::createEnd(), and SimulationAreaIterator::regionOfInterest.
Referenced by iterateOverRegionOfInterest().
|
protected |
Calculate global index from two axis indices.
Definition at line 377 of file IDetector.cpp.
References axis(), rank(), IAxis::size(), and totalSize().
Referenced by RectangularDetector::indexOfSpecular(), and SphericalDetector::indexOfSpecular().
bool IDetector::hasExplicitRegionOfInterest | ( | ) | const |
True if a region of interest is explicitly set.
Definition at line 131 of file IDetector.cpp.
References OwningVector< T >::empty(), m_axes, m_explicitROI, and OwningVector< T >::size().
|
pure virtual |
Returns index of pixel that contains the specular wavevector. If no pixel contains this specular wavevector, the number of pixels is returned. This corresponds to an overflow index.
Implemented in SphericalDetector, and RectangularDetector.
Referenced by ISimulation2D::generateElements().
void IDetector::iterateOverNonMaskedPoints | ( | std::function< void(const_iterator)> | func | ) | const |
Iterate over all non-masked points within "region of interest". If no region of interest is explicitly defined, then the whole detector is taken as "region of interest".
Definition at line 252 of file IDetector.cpp.
References beginNonMaskedPoints(), endNonMaskedPoints(), and rank().
Referenced by active_indices(), applyDetectorResolution(), and numberOfElements().
void IDetector::iterateOverRegionOfInterest | ( | std::function< void(const_iterator)> | func | ) | const |
Iterate over all points within "region of interest", no matter whether they are masked or not. If no region of interest is explicitly defined, then the whole detector is taken as "region of interest".
Definition at line 243 of file IDetector.cpp.
References beginRegionOfInterestPoints(), endRegionOfInterestPoints(), and rank().
Referenced by ScatteringSimulation::intensityMapSize().
void IDetector::maskAll | ( | ) |
Put the mask for all detector channels (i.e. exclude whole detector from the analysis)
Definition at line 365 of file IDetector.cpp.
References addMask(), and rank().
Referenced by ISimulation2D::maskAll().
|
overridevirtual |
Returns all children.
Reimplemented from INode.
Definition at line 238 of file IDetector.cpp.
References m_detector_resolution, and m_polAnalyzer.
|
inherited |
Returns all descendants.
Definition at line 61 of file INode.cpp.
References INode::nodeChildren().
size_t IDetector::numberOfElements | ( | ) | const |
Returns number of simulation elements.
Definition at line 220 of file IDetector.cpp.
References iterateOverNonMaskedPoints().
|
pure virtual |
Implemented in SphericalDetector, and RectangularDetector.
Referenced by OffspecSimulation::createCoordSystem().
|
inlinevirtualinherited |
Returns the parameter definitions, to be hard-coded in each leaf class.
Reimplemented in ConstantBackground, GaussSphere, FuzzySphere, RotationEuler, RotationZ, RotationY, RotationX, Crystal, Layer, HexagonalLattice2D, SquareLattice2D, BasicLattice2D, LayerRoughness, TruncatedSpheroid, TruncatedSphere, TruncatedCube, Spheroid, Sphere, SawtoothRippleLorentz, SawtoothRippleGauss, SawtoothRippleBox, Pyramid6, Pyramid4, Pyramid3, Pyramid2, Prism6, Prism3, PlatonicTetrahedron, PlatonicOctahedron, LongBoxLorentz, LongBoxGauss, Icosahedron, HorizontalCylinder, HollowSphere, HemiEllipsoid, EllipsoidalCylinder, Dodecahedron, Cylinder, CosineRippleLorentz, CosineRippleGauss, CosineRippleBox, Cone, CantellatedCube, Box, Bipyramid4, BarLorentz, BarGauss, Profile2DVoigt, Profile2DCone, Profile2DGate, Profile2DGauss, Profile2DCauchy, Profile1DVoigt, Profile1DCosine, Profile1DTriangle, Profile1DGate, Profile1DGauss, Profile1DCauchy, MisesGaussPeakShape, MisesFisherGaussPeakShape, LorentzFisherPeakShape, GaussFisherPeakShape, IsotropicLorentzPeakShape, IsotropicGaussPeakShape, ParticleLayout, InterferenceTwin, InterferenceRadialParaCrystal, InterferenceHardDisk, Interference2DSuperLattice, Interference2DParaCrystal, Interference1DLattice, DistributionTrapezoid, DistributionCosine, DistributionLogNormal, DistributionGaussian, DistributionLorentz, DistributionGate, ResolutionFunction2DGaussian, PolFilter, FootprintSquare, and FootprintGauss.
Definition at line 51 of file INode.h.
Referenced by INode::checkNodeArgs(), and IFormFactor::pythonConstructor().
size_t IDetector::rank | ( | ) | const |
Returns number of defined axes.
Definition at line 64 of file IDetector.cpp.
References m_axes, and OwningVector< T >::size().
Referenced by addDetAxis(), axis(), axisBinIndex(), boundsOfExplicitRegionOfInterest(), createDetectorMap(), getGlobalIndex(), RectangularDetector::indexOfSpecular(), SphericalDetector::indexOfSpecular(), iterateOverNonMaskedPoints(), iterateOverRegionOfInterest(), maskAll(), setRegionOfInterest(), totalSize(), and OffspecSimulation::transferDetectorImage().
std::pair< double, double > IDetector::regionOfInterestBounds | ( | size_t | iAxis | ) | const |
The lower and upper bound of the region of interest. If no region of interest is explicitly defined, then the whole detector is taken as "region of interest".
Definition at line 227 of file IDetector.cpp.
References ASSERT, boundsOfExplicitRegionOfInterest(), m_axes, and OwningVector< T >::size().
Referenced by axesClippedToRegionOfInterest(), createDetectorMap(), and RectangularDetector::regionOfInterestPixel().
size_t IDetector::regionOfInterestIndexToDetectorIndex | ( | size_t | regionOfInterestIndex | ) | const |
Convert an index of the region of interest to an index of the detector. If no region of interest is set, then the index stays unmodified (since ROI == detector area).
Definition at line 281 of file IDetector.cpp.
References m_explicitROI.
Referenced by SimulationAreaIterator::detectorIndex(), and SimulationAreaIterator::isMasked().
void IDetector::resetRegionOfInterest | ( | ) |
Resets region of interest making whole detector plane available for the simulation.
Definition at line 164 of file IDetector.cpp.
References m_explicitROI.
|
pure virtual |
Implemented in SphericalDetector, and RectangularDetector.
void IDetector::setAnalyzer | ( | R3 | direction, |
double | efficiency, | ||
double | total_transmission | ||
) |
Sets the polarization analyzer characteristics of the detector.
Definition at line 147 of file IDetector.cpp.
References m_polAnalyzer.
|
inlinevirtual |
Inits detector with the beam settings.
Reimplemented in RectangularDetector.
Definition at line 67 of file IDetector.h.
void IDetector::setDetectorParameters | ( | size_t | n_x, |
double | x_min, | ||
double | x_max, | ||
size_t | n_y, | ||
double | y_min, | ||
double | y_max | ||
) |
Sets equidistant axes.
Definition at line 329 of file IDetector.cpp.
References addDetAxis(), axisName(), and clear().
Referenced by RectangularDetector::RectangularDetector(), and SphericalDetector::SphericalDetector().
void IDetector::setDetectorResolution | ( | const IDetectorResolution & | p_detector_resolution | ) |
Sets the detector resolution.
Definition at line 152 of file IDetector.cpp.
References IDetectorResolution::clone(), and m_detector_resolution.
Referenced by IDetector(), and setResolutionFunction().
void IDetector::setRegionOfInterest | ( | double | xlow, |
double | ylow, | ||
double | xup, | ||
double | yup | ||
) |
Sets rectangular region of interest with lower left and upper right corners defined.
Definition at line 337 of file IDetector.cpp.
References ASSERT, axis(), m_explicitROI, and rank().
Referenced by ISimulation2D::setRegionOfInterest().
void IDetector::setResolutionFunction | ( | const IResolutionFunction2D & | resFunc | ) |
Definition at line 158 of file IDetector.cpp.
References setDetectorResolution().
|
protectedvirtual |
Return 0 if no ROI has been explicitly set. Size means number of data points.
Definition at line 94 of file IDetector.cpp.
References m_axes, m_explicitROI, and OwningVector< T >::size().
Referenced by sizeOfRegionOfInterest().
size_t IDetector::sizeOfRegionOfInterest | ( | ) | const |
The size of the "Region of Interest". Same as totalSize() if no region of interest has been explicitly set.
Definition at line 125 of file IDetector.cpp.
References sizeOfExplicitRegionOfInterest(), and totalSize().
Referenced by SimulationAreaIterator::createEnd().
size_t IDetector::totalSize | ( | ) | const |
Returns total number of pixels. Any region of interest is not taken into account.
Definition at line 114 of file IDetector.cpp.
References m_axes, and rank().
Referenced by getGlobalIndex(), RectangularDetector::indexOfSpecular(), SphericalDetector::indexOfSpecular(), and sizeOfRegionOfInterest().
|
inlinevirtualinherited |
Used for Python overriding of clone (see swig/tweaks.py)
Definition at line 32 of file ICloneable.h.
|
private |
Definition at line 248 of file IDetector.h.
Referenced by addDetAxis(), axesClippedToRegionOfInterest(), axis(), axisBinIndex(), clear(), hasExplicitRegionOfInterest(), rank(), regionOfInterestBounds(), sizeOfExplicitRegionOfInterest(), and totalSize().
|
private |
Definition at line 251 of file IDetector.h.
Referenced by addDetAxis(), addMask(), and detectorMask().
|
private |
Definition at line 250 of file IDetector.h.
Referenced by IDetector(), applyDetectorResolution(), createDetectorIntensity(), detectorResolution(), nodeChildren(), and setDetectorResolution().
|
protected |
an explicitly defined region of interest. Empty if no ROI has been defined. Vector index corresponds to axis index in m_axes
Definition at line 242 of file IDetector.h.
Referenced by boundsOfExplicitRegionOfInterest(), detectorIndexToRegionOfInterestIndex(), hasExplicitRegionOfInterest(), regionOfInterestIndexToDetectorIndex(), resetRegionOfInterest(), setRegionOfInterest(), and sizeOfExplicitRegionOfInterest().
|
protectedinherited |
Definition at line 63 of file INode.h.
Referenced by IFootprintFactor::IFootprintFactor(), INode::checkNodeArgs(), IProfile1D::pythonConstructor(), IProfile2D::pythonConstructor(), IFormFactor::pythonConstructor(), Profile1DVoigt::pythonConstructor(), and Profile2DVoigt::pythonConstructor().
|
private |
Definition at line 249 of file IDetector.h.
Referenced by analyzer(), nodeChildren(), and setAnalyzer().