BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
RectangularDetector Class Reference
Inheritance diagram for RectangularDetector:
Collaboration diagram for RectangularDetector:

Public Types

enum  EDetectorArrangement {
  GENERIC , PERPENDICULAR_TO_SAMPLE , PERPENDICULAR_TO_DIRECT_BEAM , PERPENDICULAR_TO_REFLECTED_BEAM ,
  PERPENDICULAR_TO_REFLECTED_BEAM_DPOS
}
 
using const_iterator = const SimulationAreaIterator &
 

Public Member Functions

 RectangularDetector (size_t nxbins, double width, size_t nybins, double height)
 
 RectangularDetector (const RectangularDetector &other)
 
RectangularDetectorclone () const override
 
void accept (INodeVisitor *visitor) const final
 
 ~RectangularDetector ()
 
void init (const Beam &beam) override
 
void setPosition (const kvector_t normal_to_detector, double u0, double v0, const kvector_t direction=kvector_t(0.0, -1.0, 0.0))
 
void setPerpendicularToSampleX (double distance, double u0, double v0)
 
void setPerpendicularToDirectBeam (double distance, double u0, double v0)
 
void setPerpendicularToReflectedBeam (double distance, double u0=0.0, double v0=0.0)
 
void setDirectBeamPosition (double u0, double v0)
 
double getWidth () const
 
double getHeight () const
 
size_t getNbinsX () const
 
size_t getNbinsY () const
 
kvector_t getNormalVector () const
 
double getU0 () const
 
double getV0 () const
 
kvector_t getDirectionVector () const
 
double getDistance () const
 
double getDirectBeamU0 () const
 
double getDirectBeamV0 () const
 
EDetectorArrangement getDetectorArrangment () const
 
Axes::Units defaultAxesUnits () const override
 
RectangularPixelregionOfInterestPixel () const
 
void setDetectorParameters (size_t n_x, double x_min, double x_max, size_t n_y, double y_min, double y_max)
 
void removeMasks ()
 
const DetectorMaskdetectorMask () const override
 
void addMask (const IShape2D &shape, bool mask_value=true)
 
void maskAll ()
 
const RegionOfInterestregionOfInterest () const override
 
void setRegionOfInterest (double xlow, double ylow, double xup, double yup)
 
void resetRegionOfInterest () override
 
std::vector< size_t > active_indices () const
 
std::unique_ptr< DetectorContextcreateContext () const
 
void addAxis (const IAxis &axis)
 
const IAxisgetAxis (size_t index) const
 
size_t dimension () const
 
size_t axisBinIndex (size_t index, size_t selected_axis) const
 
size_t totalSize () const
 
void setAnalyzerProperties (const kvector_t direction, double efficiency, double total_transmission)
 
void setDetectorResolution (const IDetectorResolution &p_detector_resolution)
 
void setResolutionFunction (const IResolutionFunction2D &resFunc)
 
void applyDetectorResolution (OutputData< double > *p_intensity_map) const
 
void removeDetectorResolution ()
 
const IDetectorResolutiondetectorResolution () const
 
std::unique_ptr< OutputData< double > > createDetectorMap () const
 
const DetectionPropertiesdetectionProperties () const
 
OutputData< double > * createDetectorIntensity (const std::vector< SimulationElement > &elements) const
 
size_t numberOfSimulationElements () const
 
std::vector< const INode * > getChildren () const override
 
void iterate (std::function< void(const_iterator)> func, bool visit_masks=false) const
 
virtual void transferToCPP ()
 
virtual std::string treeToString () const
 
void registerChild (INode *node)
 
virtual void setParent (const INode *newParent)
 
const INodeparent () const
 
INodeparent ()
 
int copyNumber (const INode *node) const
 
std::string displayName () const
 
ParameterPoolcreateParameterTree () const
 
ParameterPoolparameterPool () const
 
std::string parametersToString () const
 
RealParameterregisterParameter (const std::string &name, double *parpointer)
 
void registerVector (const std::string &base_name, kvector_t *p_vec, const std::string &units="nm")
 
void setParameterValue (const std::string &name, double value)
 
void setVectorValue (const std::string &base_name, kvector_t value)
 
RealParameterparameter (const std::string &name) const
 
virtual void onChange ()
 
void removeParameter (const std::string &name)
 
void removeVector (const std::string &base_name)
 
void setName (const std::string &name)
 
const std::string & getName () const
 

Static Public Member Functions

static std::string XComponentName (const std::string &base_name)
 
static std::string YComponentName (const std::string &base_name)
 
static std::string ZComponentName (const std::string &base_name)
 

Protected Member Functions

IPixelcreatePixel (size_t index) const override
 
std::string axisName (size_t index) const override
 
size_t indexOfSpecular (const Beam &beam) const override
 
void swapContent (RectangularDetector &other)
 
size_t getGlobalIndex (size_t x, size_t y) const
 
void clear ()
 
virtual std::unique_ptr< IAxiscreateAxis (size_t index, size_t n_bins, double min, double max) const
 

Protected Attributes

const size_t m_NP
 
std::vector< double > m_P
 

Private Member Functions

void setDistanceAndOffset (double distance, double u0, double v0)
 
void initNormalVector (const kvector_t central_k)
 
void initUandV (double alpha_i)
 
void setDataToDetectorMap (OutputData< double > &detectorMap, const std::vector< SimulationElement > &elements) const
 

Private Attributes

kvector_t m_normal_to_detector
 
double m_u0
 
double m_v0
 
kvector_t m_direction
 
double m_distance
 
double m_dbeam_u0
 
double m_dbeam_v0
 
EDetectorArrangement m_detector_arrangement
 
kvector_t m_u_unit
 
kvector_t m_v_unit
 
DetectorMask m_detector_mask
 
std::unique_ptr< RegionOfInterestm_region_of_interest
 
CloneableVector< IAxism_axes
 
DetectionProperties m_detection_properties
 
std::unique_ptr< IDetectorResolutionmP_detector_resolution
 
const INodem_parent {nullptr}
 
std::string m_name
 
std::unique_ptr< ParameterPoolm_pool
 

Detailed Description

A flat rectangular detector with axes and resolution function.

Definition at line 26 of file RectangularDetector.h.

Member Typedef Documentation

◆ const_iterator

Definition at line 38 of file IDetector.h.

Member Enumeration Documentation

◆ EDetectorArrangement

Enumerator
GENERIC 
PERPENDICULAR_TO_SAMPLE 
PERPENDICULAR_TO_DIRECT_BEAM 
PERPENDICULAR_TO_REFLECTED_BEAM 
PERPENDICULAR_TO_REFLECTED_BEAM_DPOS 

Definition at line 29 of file RectangularDetector.h.

Constructor & Destructor Documentation

◆ RectangularDetector() [1/2]

RectangularDetector::RectangularDetector ( size_t  nxbins,
double  width,
size_t  nybins,
double  height 
)

Rectangular detector constructor.

Parameters
nxbinsNumber of bins (pixels) in x-direction
widthWidth of the detector in mm along x-direction
nybinsNumber of bins (pixels) in y-direction
heightHeight of the detector in mm along y-direction

Definition at line 24 of file RectangularDetector.cpp.

25  : m_u0(0.0), m_v0(0.0), m_direction(kvector_t(0.0, -1.0, 0.0)), m_distance(0.0),
27 {
28  setDetectorParameters(nxbins, 0.0, width, nybins, 0.0, height);
29  setName("RectangularDetector");
30 }
BasicVector3D< double > kvector_t
Definition: Vectors3D.h:21
void setDetectorParameters(size_t n_x, double x_min, double x_max, size_t n_y, double y_min, double y_max)
Sets detector parameters using angle ranges.
Definition: IDetector2D.cpp:35
void setName(const std::string &name)
double m_dbeam_v0
position of direct beam in detector coordinates
double m_v0
position of normal vector hitting point in detector coordinates
kvector_t m_direction
direction vector of detector coordinate system
double m_distance
distance from sample origin to the detector plane
EDetectorArrangement m_detector_arrangement

References anonymous_namespace{BoxCompositionBuilder.cpp}::height, IDetector2D::setDetectorParameters(), IParameterized::setName(), and anonymous_namespace{BoxCompositionBuilder.cpp}::width.

Referenced by clone().

Here is the call graph for this function:

◆ RectangularDetector() [2/2]

RectangularDetector::RectangularDetector ( const RectangularDetector other)

Definition at line 32 of file RectangularDetector.cpp.

34  m_v0(other.m_v0), m_direction(other.m_direction), m_distance(other.m_distance),
37  m_v_unit(other.m_v_unit)
38 {
39  setName("RectangularDetector");
40 }

References IParameterized::setName().

Here is the call graph for this function:

◆ ~RectangularDetector()

RectangularDetector::~RectangularDetector ( )
default

Member Function Documentation

◆ clone()

RectangularDetector * RectangularDetector::clone ( ) const
overridevirtual

Implements IDetector2D.

Definition at line 44 of file RectangularDetector.cpp.

45 {
46  return new RectangularDetector(*this);
47 }
RectangularDetector(size_t nxbins, double width, size_t nybins, double height)
Rectangular detector constructor.

References RectangularDetector().

Here is the call graph for this function:

◆ accept()

void RectangularDetector::accept ( INodeVisitor visitor) const
inlinefinalvirtual

Calls the INodeVisitor's visit method.

Implements INode.

Definition at line 48 of file RectangularDetector.h.

48 { visitor->visit(this); }
virtual void visit(const BasicLattice *)
Definition: INodeVisitor.h:154

◆ init()

void RectangularDetector::init ( const Beam )
overridevirtual

Inits detector with the beam settings.

Reimplemented from IDetector.

Definition at line 49 of file RectangularDetector.cpp.

50 {
51  double alpha_i = beam.getAlpha();
52  kvector_t central_k = beam.getCentralK();
53  initNormalVector(central_k);
54  initUandV(alpha_i);
55 }
void initUandV(double alpha_i)
void initNormalVector(const kvector_t central_k)

References Beam::getAlpha(), Beam::getCentralK(), initNormalVector(), and initUandV().

Here is the call graph for this function:

◆ setPosition()

void RectangularDetector::setPosition ( const kvector_t  normal_to_detector,
double  u0,
double  v0,
const kvector_t  direction = kvector_t(0.0, -1.0, 0.0) 
)

Definition at line 57 of file RectangularDetector.cpp.

59 {
61  m_normal_to_detector = normal_to_detector;
63  m_u0 = u0;
64  m_v0 = v0;
65  m_direction = direction;
66 }
double mag() const
Returns magnitude of the vector.

References GENERIC, m_detector_arrangement, m_direction, m_distance, m_normal_to_detector, m_u0, m_v0, and BasicVector3D< T >::mag().

Referenced by StandardSimulations::RectDetectorGeneric().

Here is the call graph for this function:

◆ setPerpendicularToSampleX()

void RectangularDetector::setPerpendicularToSampleX ( double  distance,
double  u0,
double  v0 
)

Definition at line 68 of file RectangularDetector.cpp.

69 {
71  setDistanceAndOffset(distance, u0, v0);
72 }
void setDistanceAndOffset(double distance, double u0, double v0)

References m_detector_arrangement, PERPENDICULAR_TO_SAMPLE, and setDistanceAndOffset().

Referenced by StandardSimulations::RectDetectorPerpToSample().

Here is the call graph for this function:

◆ setPerpendicularToDirectBeam()

void RectangularDetector::setPerpendicularToDirectBeam ( double  distance,
double  u0,
double  v0 
)

Definition at line 74 of file RectangularDetector.cpp.

75 {
77  setDistanceAndOffset(distance, u0, v0);
78 }

References m_detector_arrangement, PERPENDICULAR_TO_DIRECT_BEAM, and setDistanceAndOffset().

Referenced by StandardSimulations::RectDetectorPerpToDirectBeam().

Here is the call graph for this function:

◆ setPerpendicularToReflectedBeam()

void RectangularDetector::setPerpendicularToReflectedBeam ( double  distance,
double  u0 = 0.0,
double  v0 = 0.0 
)

Definition at line 80 of file RectangularDetector.cpp.

References m_detector_arrangement, PERPENDICULAR_TO_REFLECTED_BEAM, and setDistanceAndOffset().

Referenced by StandardSimulations::RectDetectorPerpToReflectedBeam(), and StandardSimulations::RectDetectorPerpToReflectedBeamDpos().

Here is the call graph for this function:

◆ setDirectBeamPosition()

void RectangularDetector::setDirectBeamPosition ( double  u0,
double  v0 
)

◆ getWidth()

double RectangularDetector::getWidth ( ) const

Definition at line 93 of file RectangularDetector.cpp.

94 {
95  const IAxis& axis = getAxis(0);
96  return axis.getMax() - axis.getMin();
97 }
Interface for one-dimensional axes.
Definition: IAxis.h:25
virtual double getMin() const =0
Returns value of first point of axis.
virtual double getMax() const =0
Returns value of last point of axis.
const IAxis & getAxis(size_t index) const
Definition: IDetector.cpp:54

References IDetector::getAxis(), IAxis::getMax(), and IAxis::getMin().

Referenced by regionOfInterestPixel().

Here is the call graph for this function:

◆ getHeight()

double RectangularDetector::getHeight ( ) const

Definition at line 99 of file RectangularDetector.cpp.

100 {
101  const IAxis& axis = getAxis(1);
102  return axis.getMax() - axis.getMin();
103 }

References IDetector::getAxis(), IAxis::getMax(), and IAxis::getMin().

Referenced by regionOfInterestPixel().

Here is the call graph for this function:

◆ getNbinsX()

size_t RectangularDetector::getNbinsX ( ) const

Definition at line 105 of file RectangularDetector.cpp.

106 {
107  return getAxis(0).size();
108 }
virtual size_t size() const =0
retrieve the number of bins

References IDetector::getAxis(), and IAxis::size().

Here is the call graph for this function:

◆ getNbinsY()

size_t RectangularDetector::getNbinsY ( ) const

Definition at line 110 of file RectangularDetector.cpp.

111 {
112  return getAxis(1).size();
113 }

References IDetector::getAxis(), and IAxis::size().

Here is the call graph for this function:

◆ getNormalVector()

kvector_t RectangularDetector::getNormalVector ( ) const

Definition at line 115 of file RectangularDetector.cpp.

116 {
117  return m_normal_to_detector;
118 }

References m_normal_to_detector.

◆ getU0()

double RectangularDetector::getU0 ( ) const

Definition at line 120 of file RectangularDetector.cpp.

121 {
122  return m_u0;
123 }

References m_u0.

◆ getV0()

double RectangularDetector::getV0 ( ) const

Definition at line 125 of file RectangularDetector.cpp.

126 {
127  return m_v0;
128 }

References m_v0.

◆ getDirectionVector()

kvector_t RectangularDetector::getDirectionVector ( ) const

Definition at line 130 of file RectangularDetector.cpp.

131 {
132  return m_direction;
133 }

References m_direction.

◆ getDistance()

double RectangularDetector::getDistance ( ) const

Definition at line 135 of file RectangularDetector.cpp.

136 {
137  return m_distance;
138 }

References m_distance.

◆ getDirectBeamU0()

double RectangularDetector::getDirectBeamU0 ( ) const

Definition at line 140 of file RectangularDetector.cpp.

141 {
142  return m_dbeam_u0;
143 }

References m_dbeam_u0.

◆ getDirectBeamV0()

double RectangularDetector::getDirectBeamV0 ( ) const

Definition at line 145 of file RectangularDetector.cpp.

146 {
147  return m_dbeam_v0;
148 }

References m_dbeam_v0.

◆ getDetectorArrangment()

RectangularDetector::EDetectorArrangement RectangularDetector::getDetectorArrangment ( ) const

Definition at line 150 of file RectangularDetector.cpp.

151 {
152  return m_detector_arrangement;
153 }

References m_detector_arrangement.

◆ defaultAxesUnits()

Axes::Units RectangularDetector::defaultAxesUnits ( ) const
overridevirtual

return default axes units

Reimplemented from IDetector.

Definition at line 155 of file RectangularDetector.cpp.

156 {
157  return Axes::Units::MM;
158 }

◆ regionOfInterestPixel()

RectangularPixel * RectangularDetector::regionOfInterestPixel ( ) const

Definition at line 160 of file RectangularDetector.cpp.

161 {
162  const IAxis& u_axis = getAxis(0);
163  const IAxis& v_axis = getAxis(1);
164  double u_min, v_min, width, height;
165  auto p_roi = regionOfInterest();
166  if (p_roi) {
167  u_min = p_roi->getXlow();
168  v_min = p_roi->getYlow();
169  width = p_roi->getXup() - p_roi->getXlow();
170  height = p_roi->getYup() - p_roi->getYlow();
171  } else {
172  u_min = u_axis.getMin();
173  v_min = v_axis.getMin();
174  width = getWidth();
175  height = getHeight();
176  }
177  const kvector_t corner_position(m_normal_to_detector + (u_min - m_u0) * m_u_unit
178  + (v_min - m_v0) * m_v_unit);
179  const kvector_t uaxis_vector = width * m_u_unit;
180  const kvector_t vaxis_vector = height * m_v_unit;
181  return new RectangularPixel(corner_position, uaxis_vector, vaxis_vector);
182 }
const RegionOfInterest * regionOfInterest() const override
Returns region of interest if exists.
Definition: IDetector2D.cpp:43
A pixel in a RectangularDetector.

References IDetector::getAxis(), getHeight(), IAxis::getMin(), getWidth(), anonymous_namespace{BoxCompositionBuilder.cpp}::height, m_normal_to_detector, m_u0, m_u_unit, m_v0, m_v_unit, IDetector2D::regionOfInterest(), and anonymous_namespace{BoxCompositionBuilder.cpp}::width.

Referenced by RectangularConverter::RectangularConverter().

Here is the call graph for this function:

◆ createPixel()

IPixel * RectangularDetector::createPixel ( size_t  index) const
overrideprotectedvirtual

Creates an IPixel for the given OutputData object and index.

Implements IDetector2D.

Definition at line 184 of file RectangularDetector.cpp.

185 {
186  const IAxis& u_axis = getAxis(0);
187  const IAxis& v_axis = getAxis(1);
188  const size_t u_index = axisBinIndex(index, 0);
189  const size_t v_index = axisBinIndex(index, 1);
190 
191  const Bin1D u_bin = u_axis.getBin(u_index);
192  const Bin1D v_bin = v_axis.getBin(v_index);
193  const kvector_t corner_position(m_normal_to_detector + (u_bin.m_lower - m_u0) * m_u_unit
194  + (v_bin.m_lower - m_v0) * m_v_unit);
195  const kvector_t width = u_bin.getBinSize() * m_u_unit;
196  const kvector_t height = v_bin.getBinSize() * m_v_unit;
197  return new RectangularPixel(corner_position, width, height);
198 }
virtual Bin1D getBin(size_t index) const =0
retrieve a 1d bin for the given index
size_t axisBinIndex(size_t index, size_t selected_axis) const
Calculate axis index for given global index.
Definition: IDetector.cpp:61
Definition: Bin.h:20
double getBinSize() const
Definition: Bin.h:26
double m_lower
lower bound of the bin
Definition: Bin.h:23

References IDetector::axisBinIndex(), IDetector::getAxis(), IAxis::getBin(), Bin1D::getBinSize(), anonymous_namespace{BoxCompositionBuilder.cpp}::height, Bin1D::m_lower, m_normal_to_detector, m_u0, m_u_unit, m_v0, m_v_unit, and anonymous_namespace{BoxCompositionBuilder.cpp}::width.

Here is the call graph for this function:

◆ axisName()

std::string RectangularDetector::axisName ( size_t  index) const
overrideprotectedvirtual

Returns the name for the axis with given index.

Implements IDetector.

Definition at line 200 of file RectangularDetector.cpp.

201 {
202  switch (index) {
203  case 0:
204  return "u";
205  case 1:
206  return "v";
207  default:
209  "RectangularDetector::getAxisName(size_t index) -> Error! index > 1");
210  }
211 }

◆ indexOfSpecular()

size_t RectangularDetector::indexOfSpecular ( const Beam beam) const
overrideprotectedvirtual

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.

Implements IDetector2D.

Definition at line 213 of file RectangularDetector.cpp.

214 {
215  if (dimension() != 2)
216  return totalSize();
217  const double alpha = beam.getAlpha();
218  const double phi = beam.getPhi();
219  const kvector_t k_spec = vecOfLambdaAlphaPhi(beam.getWavelength(), alpha, phi);
220  const kvector_t normal_unit = m_normal_to_detector.unit();
221  const double kd = k_spec.dot(normal_unit);
222  if (kd <= 0.0)
223  return totalSize();
224  ASSERT(m_distance != 0);
225  const kvector_t rpix = k_spec * (m_distance / kd);
226  const double u = rpix.dot(m_u_unit) + m_u0;
227  const double v = rpix.dot(m_v_unit) + m_v0;
228  const IAxis& u_axis = getAxis(0); // the x axis, GISAS only
229  const IAxis& v_axis = getAxis(1); // the y axis, in reflectometer plane
230  if (!u_axis.contains(u) || !v_axis.contains(v))
231  return totalSize();
232  return getGlobalIndex(u_axis.findClosestIndex(u), v_axis.findClosestIndex(v));
233 }
#define ASSERT(condition)
Definition: Assert.h:26
BasicVector3D< double > vecOfLambdaAlphaPhi(double _lambda, double _alpha, double _phi)
Creates a vector<double> as a wavevector with given wavelength and angles.
auto dot(const BasicVector3D< U > &v) const
Returns dot product of vectors (antilinear in the first [=self] argument).
BasicVector3D< T > unit() const
Returns unit vector in direction of this. Throws for null vector.
double getWavelength() const
Definition: Beam.h:69
double getAlpha() const
Definition: Beam.h:70
double getPhi() const
Definition: Beam.h:71
virtual bool contains(double value) const
Returns true if axis contains given point.
Definition: IAxis.cpp:40
virtual size_t findClosestIndex(double value) const =0
find bin index which is best match for given value
size_t getGlobalIndex(size_t x, size_t y) const
Calculate global index from two axis indices.
Definition: IDetector2D.cpp:98
size_t dimension() const
Returns actual dimensionality of the detector (number of defined axes)
Definition: IDetector.cpp:44
size_t totalSize() const
Returns total number of pixels.
Definition: IDetector.cpp:87

References ASSERT, IAxis::contains(), IDetector::dimension(), BasicVector3D< T >::dot(), IAxis::findClosestIndex(), Beam::getAlpha(), IDetector::getAxis(), IDetector2D::getGlobalIndex(), Beam::getPhi(), Beam::getWavelength(), m_distance, m_normal_to_detector, m_u0, m_u_unit, m_v0, m_v_unit, IDetector::totalSize(), BasicVector3D< T >::unit(), and vecOfLambdaAlphaPhi().

Here is the call graph for this function:

◆ swapContent()

void RectangularDetector::swapContent ( RectangularDetector other)
protected

swap function

◆ setDistanceAndOffset()

void RectangularDetector::setDistanceAndOffset ( double  distance,
double  u0,
double  v0 
)
private

Definition at line 235 of file RectangularDetector.cpp.

236 {
237  if (distance <= 0.0) {
238  std::ostringstream message;
239  message << "RectangularDetector::setPerpendicularToSample() -> Error. "
240  << "Distance to sample can't be negative or zero";
241  throw Exceptions::LogicErrorException(message.str());
242  }
243  m_distance = distance;
244  m_u0 = u0;
245  m_v0 = v0;
246 }

References m_distance, m_u0, and m_v0.

Referenced by setPerpendicularToDirectBeam(), setPerpendicularToReflectedBeam(), and setPerpendicularToSampleX().

◆ initNormalVector()

void RectangularDetector::initNormalVector ( const kvector_t  central_k)
private

Definition at line 248 of file RectangularDetector.cpp.

249 {
250  kvector_t central_k_unit = central_k.unit();
251 
253  // do nothing
254  }
255 
258  }
259 
261  m_normal_to_detector = m_distance * central_k_unit;
262  }
263 
266  m_normal_to_detector = m_distance * central_k_unit;
268  }
269 
270  else {
272  "RectangularDetector::init() -> Unknown detector arrangement");
273  }
274 }
T z() const
Returns z-component in cartesian coordinate system.
Definition: BasicVector3D.h:68
void setZ(const T &a)
Sets z-component in cartesian coordinate system.
Definition: BasicVector3D.h:75

References GENERIC, m_detector_arrangement, m_distance, m_normal_to_detector, PERPENDICULAR_TO_DIRECT_BEAM, PERPENDICULAR_TO_REFLECTED_BEAM, PERPENDICULAR_TO_REFLECTED_BEAM_DPOS, PERPENDICULAR_TO_SAMPLE, BasicVector3D< T >::setZ(), BasicVector3D< T >::unit(), and BasicVector3D< T >::z().

Referenced by init().

Here is the call graph for this function:

◆ initUandV()

void RectangularDetector::initUandV ( double  alpha_i)
private

Definition at line 276 of file RectangularDetector.cpp.

277 {
279  kvector_t u_direction =
281  m_u_unit = u_direction.unit();
283 
285  kvector_t z(0.0, 0.0, 1.0);
286  kvector_t normal_unit = m_normal_to_detector.unit();
287  kvector_t zp = z - z.dot(normal_unit) * normal_unit;
288  double uz = zp.dot(m_u_unit) / zp.mag();
289  double vz = zp.dot(m_v_unit) / zp.mag();
290  m_u0 = m_dbeam_u0 + m_distance * std::tan(2 * alpha_i) * uz;
291  m_v0 = m_dbeam_v0 + m_distance * std::tan(2 * alpha_i) * vz;
292  }
293 }
auto cross(const BasicVector3D< U > &v) const
Returns cross product of vectors (linear in both arguments).

References BasicVector3D< T >::cross(), BasicVector3D< T >::dot(), m_dbeam_u0, m_dbeam_v0, m_detector_arrangement, m_direction, m_distance, m_normal_to_detector, m_u0, m_u_unit, m_v0, m_v_unit, BasicVector3D< T >::mag(), PERPENDICULAR_TO_REFLECTED_BEAM_DPOS, and BasicVector3D< T >::unit().

Referenced by init().

Here is the call graph for this function:

◆ setDetectorParameters()

void IDetector2D::setDetectorParameters ( size_t  n_x,
double  x_min,
double  x_max,
size_t  n_y,
double  y_min,
double  y_max 
)
inherited

Sets detector parameters using angle ranges.

Definition at line 35 of file IDetector2D.cpp.

37 {
38  clear();
39  addAxis(*createAxis(0, n_x, x_min, x_max));
40  addAxis(*createAxis(1, n_y, y_min, y_max));
41 }
void clear()
Definition: IDetector.cpp:49
virtual std::unique_ptr< IAxis > createAxis(size_t index, size_t n_bins, double min, double max) const
Generates an axis with correct name and default binning for given index.
Definition: IDetector.cpp:76
void addAxis(const IAxis &axis)
Definition: IDetector.cpp:39

References IDetector::addAxis(), IDetector::clear(), and IDetector::createAxis().

Referenced by IsGISAXSDetector::IsGISAXSDetector(), StandardSimulations::IsGISAXSSimulation1(), StandardSimulations::IsGISAXSSimulation2(), RectangularDetector(), Simulation2D::setDetectorParameters(), and SphericalDetector::SphericalDetector().

Here is the call graph for this function:

◆ removeMasks()

void IDetector2D::removeMasks ( )
inherited

Removes all masks from the detector.

Definition at line 74 of file IDetector2D.cpp.

75 {
77 }
void removeMasks()
remove all masks and return object to initial state
DetectorMask m_detector_mask
Definition: IDetector2D.h:89

References IDetector2D::m_detector_mask, and DetectorMask::removeMasks().

Referenced by Simulation2D::removeMasks().

Here is the call graph for this function:

◆ detectorMask()

const DetectorMask * IDetector2D::detectorMask ( ) const
overridevirtualinherited

Returns detector masks container.

Implements IDetector.

Definition at line 93 of file IDetector2D.cpp.

94 {
95  return &m_detector_mask;
96 }

References IDetector2D::m_detector_mask.

◆ addMask()

void IDetector2D::addMask ( const IShape2D shape,
bool  mask_value = true 
)
inherited

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.

Parameters
shapeThe shape of mask (Rectangle, Polygon, Line, Ellipse)
mask_valueThe value of mask

Definition at line 79 of file IDetector2D.cpp.

80 {
81  m_detector_mask.addMask(shape, mask_value);
83 }
void addMask(const IShape2D &shape, bool mask_value)
Add mask to the stack of detector masks.
void initMaskData(const IDetector2D &detector)
Init the map of masks for the given detector plane.

References DetectorMask::addMask(), DetectorMask::initMaskData(), and IDetector2D::m_detector_mask.

Referenced by Simulation2D::addMask(), and IDetector2D::maskAll().

Here is the call graph for this function:

◆ maskAll()

void IDetector2D::maskAll ( )
inherited

Put the mask for all detector channels (i.e. exclude whole detector from the analysis)

Definition at line 85 of file IDetector2D.cpp.

86 {
87  if (dimension() != 2)
88  return;
90  addMask(InfinitePlane(), true);
91 }
void addMask(const IShape2D &shape, bool mask_value=true)
Adds mask of given shape to the stack of detector masks.
Definition: IDetector2D.cpp:79
The infinite plane is used for masking everything once and forever.
Definition: InfinitePlane.h:24

References IDetector2D::addMask(), IDetector::dimension(), IDetector2D::m_detector_mask, and DetectorMask::removeMasks().

Referenced by Simulation2D::maskAll().

Here is the call graph for this function:

◆ regionOfInterest()

const RegionOfInterest * IDetector2D::regionOfInterest ( ) const
overridevirtualinherited

Returns region of interest if exists.

Implements IDetector.

Definition at line 43 of file IDetector2D.cpp.

44 {
45  return m_region_of_interest.get();
46 }
std::unique_ptr< RegionOfInterest > m_region_of_interest
Definition: IDetector2D.h:90

References IDetector2D::m_region_of_interest.

Referenced by OffSpecularConverter::addDetectorYAxis(), IDetector2D::IDetector2D(), and regionOfInterestPixel().

◆ setRegionOfInterest()

void IDetector2D::setRegionOfInterest ( double  xlow,
double  ylow,
double  xup,
double  yup 
)
inherited

Sets rectangular region of interest with lower left and upper right corners defined.

Definition at line 48 of file IDetector2D.cpp.

49 {
50  m_region_of_interest = std::make_unique<RegionOfInterest>(*this, xlow, ylow, xup, yup);
52 }

References DetectorMask::initMaskData(), IDetector2D::m_detector_mask, and IDetector2D::m_region_of_interest.

Referenced by Simulation2D::setRegionOfInterest().

Here is the call graph for this function:

◆ resetRegionOfInterest()

void IDetector2D::resetRegionOfInterest ( )
overridevirtualinherited

Resets region of interest making whole detector plane available for the simulation.

Implements IDetector.

Definition at line 54 of file IDetector2D.cpp.

55 {
56  m_region_of_interest.reset();
58 }

References DetectorMask::initMaskData(), IDetector2D::m_detector_mask, and IDetector2D::m_region_of_interest.

Here is the call graph for this function:

◆ active_indices()

std::vector< size_t > IDetector2D::active_indices ( ) const
inherited

Returns vector of unmasked detector indices.

Definition at line 60 of file IDetector2D.cpp.

61 {
62  std::vector<size_t> result;
63  SimulationArea area(this);
64  for (SimulationArea::iterator it = area.begin(); it != area.end(); ++it)
65  result.push_back(it.detectorIndex());
66  return result;
67 }
An iterator for SimulationArea.
Holds iteration logic over active detector channels in the presence of masked areas and RegionOfInter...

References SimulationArea::begin(), and SimulationArea::end().

Referenced by DetectorContext::setup_context().

Here is the call graph for this function:

◆ createContext()

std::unique_ptr< DetectorContext > IDetector2D::createContext ( ) const
inherited

Definition at line 69 of file IDetector2D.cpp.

70 {
71  return std::make_unique<DetectorContext>(this);
72 }

Referenced by Simulation2D::prepareSimulation().

◆ getGlobalIndex()

size_t IDetector2D::getGlobalIndex ( size_t  x,
size_t  y 
) const
protectedinherited

Calculate global index from two axis indices.

Definition at line 98 of file IDetector2D.cpp.

99 {
100  if (dimension() != 2)
101  return totalSize();
102  return x * getAxis(1).size() + y;
103 }

References IDetector::dimension(), IDetector::getAxis(), IAxis::size(), and IDetector::totalSize().

Referenced by indexOfSpecular(), and SphericalDetector::indexOfSpecular().

Here is the call graph for this function:

◆ addAxis()

void IDetector::addAxis ( const IAxis axis)
inherited

Definition at line 39 of file IDetector.cpp.

40 {
41  m_axes.push_back(axis.clone());
42 }
void push_back(T *t)
virtual IAxis * clone() const =0
clone function
CloneableVector< IAxis > m_axes
Definition: IDetector.h:126

References IAxis::clone(), IDetector::m_axes, and CloneableVector< T >::push_back().

Referenced by IDetector2D::setDetectorParameters(), and SpecularDetector1D::SpecularDetector1D().

Here is the call graph for this function:

◆ getAxis()

const IAxis & IDetector::getAxis ( size_t  index) const
inherited

Definition at line 54 of file IDetector.cpp.

55 {
56  if (index < dimension())
57  return *m_axes[index];
58  throw std::runtime_error("Error in IDetector::getAxis: not so many axes in this detector.");
59 }

References IDetector::dimension(), and IDetector::m_axes.

Referenced by UnitConverterSimple::addDetectorAxis(), OffSpecularConverter::addDetectorYAxis(), IDetector::createDetectorMap(), createPixel(), SphericalDetector::createPixel(), anonymous_namespace{Simulation.cpp}::detHasSameDimensions(), IDetector2D::getGlobalIndex(), getHeight(), getNbinsX(), getNbinsY(), getWidth(), indexOfSpecular(), SphericalDetector::indexOfSpecular(), DetectorMask::initMaskData(), RegionOfInterest::RegionOfInterest(), and regionOfInterestPixel().

Here is the call graph for this function:

◆ dimension()

◆ axisBinIndex()

size_t IDetector::axisBinIndex ( size_t  index,
size_t  selected_axis 
) const
inherited

Calculate axis index for given global index.

Definition at line 61 of file IDetector.cpp.

62 {
63  const size_t dim = dimension();
64  size_t remainder(index);
65  size_t i_axis = dim;
66  for (size_t i = 0; i < dim; ++i) {
67  --i_axis;
68  if (selected_axis == i_axis)
69  return remainder % m_axes[i_axis]->size();
70  remainder /= m_axes[i_axis]->size();
71  }
72  throw std::runtime_error("IDetector::getAxisBinIndex() -> "
73  "Error! No axis with given number");
74 }

References IDetector::dimension(), and IDetector::m_axes.

Referenced by createPixel(), and SphericalDetector::createPixel().

Here is the call graph for this function:

◆ totalSize()

size_t IDetector::totalSize ( ) const
inherited

Returns total number of pixels.

Definition at line 87 of file IDetector.cpp.

88 {
89  const size_t dim = dimension();
90  if (dim == 0)
91  return 0;
92  size_t result = 1;
93  for (size_t i_axis = 0; i_axis < dim; ++i_axis)
94  result *= m_axes[i_axis]->size();
95  return result;
96 }

References IDetector::dimension(), and IDetector::m_axes.

Referenced by IDetector2D::getGlobalIndex(), IsGISAXSDetector::indexOfSpecular(), indexOfSpecular(), SphericalDetector::indexOfSpecular(), and SimulationArea::SimulationArea().

Here is the call graph for this function:

◆ setAnalyzerProperties()

void IDetector::setAnalyzerProperties ( const kvector_t  direction,
double  efficiency,
double  total_transmission 
)
inherited

Sets the polarization analyzer characteristics of the detector.

Definition at line 98 of file IDetector.cpp.

100 {
101  m_detection_properties.setAnalyzerProperties(direction, efficiency, total_transmission);
102 }
void setAnalyzerProperties(const kvector_t direction, double efficiency, double total_transmission)
Sets the polarization analyzer characteristics of the detector.
DetectionProperties m_detection_properties
Definition: IDetector.h:127

References IDetector::m_detection_properties, and DetectionProperties::setAnalyzerProperties().

Here is the call graph for this function:

◆ setDetectorResolution()

void IDetector::setDetectorResolution ( const IDetectorResolution p_detector_resolution)
inherited

Sets the detector resolution.

Definition at line 104 of file IDetector.cpp.

105 {
106  mP_detector_resolution.reset(p_detector_resolution.clone());
108 }
virtual IDetectorResolution * clone() const =0
std::unique_ptr< IDetectorResolution > mP_detector_resolution
Definition: IDetector.h:128
void registerChild(INode *node)
Definition: INode.cpp:58

References IDetectorResolution::clone(), IDetector::mP_detector_resolution, and INode::registerChild().

Referenced by IDetector::IDetector(), and IDetector::setResolutionFunction().

Here is the call graph for this function:

◆ setResolutionFunction()

void IDetector::setResolutionFunction ( const IResolutionFunction2D resFunc)
inherited

Definition at line 111 of file IDetector.cpp.

112 {
113  ConvolutionDetectorResolution convFunc(resFunc);
114  setDetectorResolution(convFunc);
115 }
Convolutes the intensity in 1 or 2 dimensions with a resolution function.
void setDetectorResolution(const IDetectorResolution &p_detector_resolution)
Sets the detector resolution.
Definition: IDetector.cpp:104

References IDetector::setDetectorResolution().

Here is the call graph for this function:

◆ applyDetectorResolution()

void IDetector::applyDetectorResolution ( OutputData< double > *  p_intensity_map) const
inherited

Applies the detector resolution to the given intensity maps.

Definition at line 117 of file IDetector.cpp.

118 {
119  if (!p_intensity_map)
120  throw std::runtime_error("IDetector::applyDetectorResolution() -> "
121  "Error! Null pointer to intensity map");
123  mP_detector_resolution->applyDetectorResolution(p_intensity_map);
124  if (detectorMask() && detectorMask()->hasMasks()) {
125  // sets amplitude in masked areas to zero
126  std::unique_ptr<OutputData<double>> buff(new OutputData<double>());
127  buff->copyShapeFrom(*p_intensity_map);
128 
129  iterate([&](const_iterator it) {
130  (*buff)[it.roiIndex()] = (*p_intensity_map)[it.roiIndex()];
131  });
132  p_intensity_map->setRawDataVector(buff->getRawDataVector());
133  }
134  }
135 }
virtual const DetectorMask * detectorMask() const =0
Returns detector masks container.
const SimulationAreaIterator & const_iterator
Definition: IDetector.h:38
void iterate(std::function< void(const_iterator)> func, bool visit_masks=false) const
Definition: IDetector.cpp:201
void setRawDataVector(const std::vector< T > &data_vector)
Sets new values to raw data vector.
Definition: OutputData.h:559

References IDetector::detectorMask(), IDetector::iterate(), IDetector::mP_detector_resolution, SimulationAreaIterator::roiIndex(), and OutputData< T >::setRawDataVector().

Referenced by IDetector::createDetectorIntensity().

Here is the call graph for this function:

◆ removeDetectorResolution()

void IDetector::removeDetectorResolution ( )
inherited

Removes detector resolution function.

Definition at line 137 of file IDetector.cpp.

138 {
139  mP_detector_resolution.reset();
140 }

References IDetector::mP_detector_resolution.

◆ detectorResolution()

const IDetectorResolution * IDetector::detectorResolution ( ) const
inherited

Returns a pointer to detector resolution object.

Definition at line 142 of file IDetector.cpp.

143 {
144  return mP_detector_resolution.get();
145 }

References IDetector::mP_detector_resolution.

Referenced by SimulationToPython::defineDetectorResolutionFunction().

◆ createDetectorMap()

std::unique_ptr< OutputData< double > > IDetector::createDetectorMap ( ) const
inherited

Returns empty detector map in given axes units.

Definition at line 162 of file IDetector.cpp.

163 {
164  const size_t dim = dimension();
165  if (dim == 0)
166  throw std::runtime_error(
167  "Error in IDetector::createDetectorMap: dimensions of the detector are undefined");
168 
169  std::unique_ptr<OutputData<double>> result(new OutputData<double>);
170  for (size_t i = 0; i < dim; ++i)
171  if (auto roi = regionOfInterest())
172  result->addAxis(*roi->clipAxisToRoi(i, getAxis(i)));
173  else
174  result->addAxis(getAxis(i));
175 
176  return result;
177 }
virtual const RegionOfInterest * regionOfInterest() const =0
Returns region of interest if exists.

References IDetector::dimension(), IDetector::getAxis(), and IDetector::regionOfInterest().

Referenced by IDetector::createDetectorIntensity().

Here is the call graph for this function:

◆ detectionProperties()

◆ createDetectorIntensity()

OutputData< double > * IDetector::createDetectorIntensity ( const std::vector< SimulationElement > &  elements) const
inherited

Returns new intensity map with detector resolution applied.

Map will be cropped to ROI if ROI is present.

Definition at line 148 of file IDetector.cpp.

149 {
150  std::unique_ptr<OutputData<double>> detectorMap(createDetectorMap());
151  if (!detectorMap)
152  throw Exceptions::RuntimeErrorException("Instrument::createDetectorIntensity:"
153  "can't create detector map.");
154 
155  setDataToDetectorMap(*detectorMap, elements);
157  applyDetectorResolution(detectorMap.get());
158 
159  return detectorMap.release();
160 }
void setDataToDetectorMap(OutputData< double > &detectorMap, const std::vector< SimulationElement > &elements) const
Definition: IDetector.cpp:179
void applyDetectorResolution(OutputData< double > *p_intensity_map) const
Applies the detector resolution to the given intensity maps.
Definition: IDetector.cpp:117
std::unique_ptr< OutputData< double > > createDetectorMap() const
Returns empty detector map in given axes units.
Definition: IDetector.cpp:162

References IDetector::applyDetectorResolution(), IDetector::createDetectorMap(), IDetector::mP_detector_resolution, and IDetector::setDataToDetectorMap().

Here is the call graph for this function:

◆ numberOfSimulationElements()

size_t IDetector::numberOfSimulationElements ( ) const
inherited

Returns number of simulation elements.

Definition at line 189 of file IDetector.cpp.

190 {
191  size_t result(0);
192  iterate([&result](const_iterator) { ++result; });
193  return result;
194 }

References IDetector::iterate().

Here is the call graph for this function:

◆ getChildren()

std::vector< const INode * > IDetector::getChildren ( ) const
overridevirtualinherited

Returns a vector of children (const).

Reimplemented from INode.

Definition at line 196 of file IDetector.cpp.

197 {
198  return std::vector<const INode*>() << &m_detection_properties << mP_detector_resolution;
199 }

References IDetector::m_detection_properties, and IDetector::mP_detector_resolution.

◆ iterate()

void IDetector::iterate ( std::function< void(const_iterator)>  func,
bool  visit_masks = false 
) const
inherited

Definition at line 201 of file IDetector.cpp.

202 {
203  if (this->dimension() == 0)
204  return;
205 
206  if (visit_masks) {
207  SimulationRoiArea area(this);
208  for (SimulationRoiArea::iterator it = area.begin(); it != area.end(); ++it)
209  func(it);
210  } else {
211  SimulationArea area(this);
212  for (SimulationArea::iterator it = area.begin(); it != area.end(); ++it)
213  func(it);
214  }
215 }
Holds iteration logic over active detector channels in the presence of ROI.

References SimulationArea::begin(), IDetector::dimension(), and SimulationArea::end().

Referenced by IDetector::applyDetectorResolution(), Simulation::convertData(), GISASSimulation::intensityMapSize(), IDetector::numberOfSimulationElements(), and IDetector::setDataToDetectorMap().

Here is the call graph for this function:

◆ clear()

void IDetector::clear ( )
protectedinherited

Definition at line 49 of file IDetector.cpp.

50 {
51  m_axes.clear();
52 }

References IDetector::m_axes.

Referenced by IDetector2D::setDetectorParameters().

◆ createAxis()

std::unique_ptr< IAxis > IDetector::createAxis ( size_t  index,
size_t  n_bins,
double  min,
double  max 
) const
protectedvirtualinherited

Generates an axis with correct name and default binning for given index.

Reimplemented in IsGISAXSDetector.

Definition at line 76 of file IDetector.cpp.

78 {
79  if (max <= min)
80  throw Exceptions::LogicErrorException("IDetector::createAxis() -> Error! max <= min");
81  if (n_bins == 0)
83  "IDetector::createAxis() -> Error! Number n_bins can't be zero.");
84  return std::make_unique<FixedBinAxis>(axisName(index), n_bins, min, max);
85 }
virtual std::string axisName(size_t index) const =0
Returns the name for the axis with given index.

References IDetector::axisName().

Referenced by IDetector2D::setDetectorParameters().

Here is the call graph for this function:

◆ setDataToDetectorMap()

void IDetector::setDataToDetectorMap ( OutputData< double > &  detectorMap,
const std::vector< SimulationElement > &  elements 
) const
privateinherited

Definition at line 179 of file IDetector.cpp.

181 {
182  if (elements.empty())
183  return;
184  iterate([&](const_iterator it) {
185  detectorMap[it.roiIndex()] = elements[it.elementIndex()].getIntensity();
186  });
187 }

References SimulationAreaIterator::elementIndex(), IDetector::iterate(), and SimulationAreaIterator::roiIndex().

Referenced by IDetector::createDetectorIntensity().

Here is the call graph for this function:

◆ transferToCPP()

virtual void ICloneable::transferToCPP ( )
inlinevirtualinherited

Used for Python overriding of clone (see swig/tweaks.py)

Definition at line 34 of file ICloneable.h.

◆ treeToString()

std::string INode::treeToString ( ) const
virtualinherited

Returns multiline string representing tree structure below the node.

Definition at line 53 of file INode.cpp.

54 {
55  return NodeUtils::nodeToString(*this);
56 }
std::string nodeToString(const INode &node)
Returns multiline string representing tree structure starting from given node.
Definition: NodeUtils.cpp:67

References NodeUtils::nodeToString().

Here is the call graph for this function:

◆ registerChild()

void INode::registerChild ( INode node)
inherited

Definition at line 58 of file INode.cpp.

59 {
60  ASSERT(node);
61  node->setParent(this);
62 }
virtual void setParent(const INode *newParent)
Definition: INode.cpp:69

References ASSERT, and INode::setParent().

Referenced by ParticleLayout::addAndRegisterAbstractParticle(), ParticleCoreShell::addAndRegisterCore(), MultiLayer::addAndRegisterInterface(), MultiLayer::addAndRegisterLayer(), ParticleCoreShell::addAndRegisterShell(), Layer::addLayout(), ParticleComposition::addParticlePointer(), Beam::Beam(), Crystal::Crystal(), IDetector::IDetector(), Simulation::initialize(), MesoCrystal::initialize(), Instrument::Instrument(), Beam::operator=(), Instrument::operator=(), Particle::Particle(), ParticleDistribution::ParticleDistribution(), IParticle::rotate(), ParticleLayout::setAndRegisterInterferenceFunction(), Simulation::setBackground(), InterferenceFunction1DLattice::setDecayFunction(), InterferenceFunction2DLattice::setDecayFunction(), Instrument::setDetector(), IDetector::setDetectorResolution(), Beam::setFootprintFactor(), Particle::setFormFactor(), InterferenceFunctionFinite3DLattice::setLattice(), InterferenceFunction2DLattice::setLattice(), InterferenceFunction2DParaCrystal::setLattice(), InterferenceFunction2DSuperLattice::setLattice(), InterferenceFunctionFinite2DLattice::setLattice(), InterferenceFunctionRadialParaCrystal::setProbabilityDistribution(), InterferenceFunction2DParaCrystal::setProbabilityDistributions(), ConvolutionDetectorResolution::setResolutionFunction(), IParticle::setRotation(), LayerInterface::setRoughness(), and InterferenceFunction2DSuperLattice::setSubstructureIFF().

Here is the call graph for this function:

◆ setParent()

void INode::setParent ( const INode newParent)
virtualinherited

Reimplemented in SampleProvider.

Definition at line 69 of file INode.cpp.

70 {
71  m_parent = newParent;
72 }
const INode * m_parent
Definition: INode.h:81

References INode::m_parent.

Referenced by INode::registerChild(), SampleProvider::setBuilder(), and SampleProvider::setParent().

◆ parent() [1/2]

const INode * INode::parent ( ) const
inherited

◆ parent() [2/2]

INode * INode::parent ( )
inherited

Definition at line 79 of file INode.cpp.

80 {
81  return const_cast<INode*>(m_parent);
82 }
Base class for tree-like structures containing parameterized objects.
Definition: INode.h:49

References INode::m_parent.

◆ copyNumber()

int INode::copyNumber ( const INode node) const
inherited

Returns copyNumber of child, which takes into account existence of children with same name.

Definition at line 84 of file INode.cpp.

85 {
86  if (node->parent() != this)
87  return -1;
88 
89  int result(-1), count(0);
90  for (auto child : getChildren()) {
91 
92  if (child == nullptr)
93  throw std::runtime_error("INode::copyNumber() -> Error. Nullptr as child.");
94 
95  if (child == node)
96  result = count;
97 
98  if (child->getName() == node->getName())
99  ++count;
100  }
101 
102  return count > 1 ? result : -1;
103 }
const INode * parent() const
Definition: INode.cpp:74
virtual std::vector< const INode * > getChildren() const
Returns a vector of children (const).
Definition: INode.cpp:64
const std::string & getName() const

References INode::getChildren(), IParameterized::getName(), and INode::parent().

Referenced by INode::displayName().

Here is the call graph for this function:

◆ displayName()

std::string INode::displayName ( ) const
inherited

Returns display name, composed from the name of node and it's copy number.

Definition at line 105 of file INode.cpp.

106 {
107  std::string result = getName();
108  if (m_parent) {
109  int index = m_parent->copyNumber(this);
110  if (index >= 0)
111  result = result + std::to_string(index);
112  }
113  return result;
114 }
int copyNumber(const INode *node) const
Returns copyNumber of child, which takes into account existence of children with same name.
Definition: INode.cpp:84

References INode::copyNumber(), IParameterized::getName(), and INode::m_parent.

Referenced by NodeUtils::nodePath(), and anonymous_namespace{NodeUtils.cpp}::nodeString().

Here is the call graph for this function:

◆ createParameterTree()

ParameterPool * INode::createParameterTree ( ) const
virtualinherited

Creates new parameter pool, with all local parameters and those of its children.

Reimplemented from IParameterized.

Definition at line 116 of file INode.cpp.

117 {
118  std::unique_ptr<ParameterPool> result(new ParameterPool);
119 
121  it.first();
122  while (!it.isDone()) {
123  const INode* child = it.getCurrent();
124  const std::string path = NodeUtils::nodePath(*child, this->parent()) + "/";
125  child->parameterPool()->copyToExternalPool(path, result.get());
126  it.next();
127  }
128 
129  return result.release();
130 }
ParameterPool * parameterPool() const
Returns pointer to the parameter pool.
Iterator through INode tree of objects.
Definition: NodeIterator.h:90
Container with parameters for IParameterized object.
Definition: ParameterPool.h:30
void copyToExternalPool(const std::string &prefix, ParameterPool *other_pool) const
Copies parameters of given pool to other pool, prepeding prefix to the parameter names.
std::string nodePath(const INode &node, const INode *root=nullptr)
Returns path composed of node's displayName, with respect to root node.
Definition: NodeUtils.cpp:82

References ParameterPool::copyToExternalPool(), NodeIterator< Strategy >::first(), NodeIterator< Strategy >::getCurrent(), NodeIterator< Strategy >::isDone(), NodeIterator< Strategy >::next(), NodeUtils::nodePath(), IParameterized::parameterPool(), and INode::parent().

Referenced by ParticleDistribution::generateParticles(), Simulation::runSimulation(), DepthProbeSimulation::validateParametrization(), OffSpecSimulation::validateParametrization(), and SpecularSimulation::validateParametrization().

Here is the call graph for this function:

◆ parameterPool()

ParameterPool* IParameterized::parameterPool ( ) const
inlineinherited

Returns pointer to the parameter pool.

Definition at line 38 of file IParameterized.h.

38 { return m_pool.get(); } // has non-const usages!
std::unique_ptr< ParameterPool > m_pool
parameter pool (kind of pointer-to-implementation)

References IParameterized::m_pool.

Referenced by pyfmt2::argumentList(), SampleBuilderNode::borrow_builder_parameters(), INode::createParameterTree(), INode::INode(), IParameterized::IParameterized(), anonymous_namespace{NodeUtils.cpp}::poolToString(), SampleBuilderNode::reset(), and IDistribution1D::setUnits().

◆ parametersToString()

std::string IParameterized::parametersToString ( ) const
inherited

Returns multiline string representing available parameters.

Definition at line 40 of file IParameterized.cpp.

41 {
42  std::ostringstream result;
43  std::unique_ptr<ParameterPool> P_pool(createParameterTree());
44  result << *P_pool << "\n";
45  return result.str();
46 }
virtual ParameterPool * createParameterTree() const
Creates new parameter pool, with all local parameters and those of its children.

References IParameterized::createParameterTree().

Here is the call graph for this function:

◆ registerParameter()

RealParameter & IParameterized::registerParameter ( const std::string &  name,
double *  parpointer 
)
inherited

Definition at line 48 of file IParameterized.cpp.

49 {
50  return m_pool->addParameter(
51  new RealParameter(name, data, getName(), [&]() -> void { onChange(); }));
52 }
virtual void onChange()
Action to be taken in inherited class when a parameter has changed.
Wraps a parameter of type double.
Definition: RealParameter.h:32

References IParameterized::getName(), IParameterized::m_pool, and IParameterized::onChange().

Referenced by BasicLattice::BasicLattice(), Beam::Beam(), CylindersInBABuilder::CylindersInBABuilder(), DetectionProperties::DetectionProperties(), HexagonalLattice::HexagonalLattice(), IInterferenceFunction::IInterferenceFunction(), INode::INode(), InterferenceFunction1DLattice::InterferenceFunction1DLattice(), InterferenceFunction2DParaCrystal::InterferenceFunction2DParaCrystal(), InterferenceFunctionHardDisk::InterferenceFunctionHardDisk(), InterferenceFunctionRadialParaCrystal::InterferenceFunctionRadialParaCrystal(), InterferenceFunctionTwin::InterferenceFunctionTwin(), Lattice2D::Lattice2D(), LayerRoughness::LayerRoughness(), MultiLayer::MultiLayer(), ParticleDistribution::ParticleDistribution(), PlainMultiLayerBySLDBuilder::PlainMultiLayerBySLDBuilder(), IParticle::registerAbundance(), ParticleLayout::registerParticleDensity(), Layer::registerThickness(), IParameterized::registerVector(), ParticleLayout::registerWeight(), ResolutionFunction2DGaussian::ResolutionFunction2DGaussian(), ResonatorBuilder::ResonatorBuilder(), Lattice2D::setRotationEnabled(), SquareLattice::SquareLattice(), and TriangularRippleBuilder::TriangularRippleBuilder().

Here is the call graph for this function:

◆ registerVector()

void IParameterized::registerVector ( const std::string &  base_name,
kvector_t p_vec,
const std::string &  units = "nm" 
)
inherited

Definition at line 54 of file IParameterized.cpp.

56 {
57  registerParameter(XComponentName(base_name), &((*p_vec)[0])).setUnit(units);
58  registerParameter(YComponentName(base_name), &((*p_vec)[1])).setUnit(units);
59  registerParameter(ZComponentName(base_name), &((*p_vec)[2])).setUnit(units);
60 }
RealParameter & registerParameter(const std::string &name, double *parpointer)
static std::string XComponentName(const std::string &base_name)
static std::string YComponentName(const std::string &base_name)
static std::string ZComponentName(const std::string &base_name)
RealParameter & setUnit(const std::string &name)

References IParameterized::registerParameter(), RealParameter::setUnit(), IParameterized::XComponentName(), IParameterized::YComponentName(), and IParameterized::ZComponentName().

Referenced by Beam::Beam(), DetectionProperties::DetectionProperties(), InterferenceFunctionTwin::InterferenceFunctionTwin(), MultiLayer::MultiLayer(), Lattice::registerBasisVectors(), and IParticle::registerPosition().

Here is the call graph for this function:

◆ setParameterValue()

void IParameterized::setParameterValue ( const std::string &  name,
double  value 
)
inherited

Definition at line 62 of file IParameterized.cpp.

63 {
64  if (name.find('*') == std::string::npos && name.find('/') == std::string::npos) {
65  m_pool->setParameterValue(name, value);
66  } else {
67  std::unique_ptr<ParameterPool> P_pool{createParameterTree()};
68  if (name.find('*') != std::string::npos)
69  P_pool->setMatchedParametersValue(name, value);
70  else
71  P_pool->setParameterValue(name, value);
72  }
73 }
int setMatchedParametersValue(const std::string &wildcards, double value)
Sets value of the nonzero parameters that match pattern ('*' allowed), or throws.

References IParameterized::createParameterTree(), IParameterized::m_pool, and ParameterPool::setMatchedParametersValue().

Referenced by AsymRippleBuilder::buildSample(), and IParameterized::setVectorValue().

Here is the call graph for this function:

◆ setVectorValue()

void IParameterized::setVectorValue ( const std::string &  base_name,
kvector_t  value 
)
inherited

Definition at line 75 of file IParameterized.cpp.

76 {
77  setParameterValue(XComponentName(base_name), value.x());
78  setParameterValue(YComponentName(base_name), value.y());
79  setParameterValue(ZComponentName(base_name), value.z());
80 }
T y() const
Returns y-component in cartesian coordinate system.
Definition: BasicVector3D.h:66
T x() const
Returns x-component in cartesian coordinate system.
Definition: BasicVector3D.h:64
void setParameterValue(const std::string &name, double value)

References IParameterized::setParameterValue(), BasicVector3D< T >::x(), IParameterized::XComponentName(), BasicVector3D< T >::y(), IParameterized::YComponentName(), BasicVector3D< T >::z(), and IParameterized::ZComponentName().

Here is the call graph for this function:

◆ parameter()

RealParameter * IParameterized::parameter ( const std::string &  name) const
inherited

◆ onChange()

◆ removeParameter()

void IParameterized::removeParameter ( const std::string &  name)
inherited

◆ removeVector()

void IParameterized::removeVector ( const std::string &  base_name)
inherited

Definition at line 93 of file IParameterized.cpp.

94 {
95  removeParameter(XComponentName(base_name));
96  removeParameter(YComponentName(base_name));
97  removeParameter(ZComponentName(base_name));
98 }
void removeParameter(const std::string &name)

References IParameterized::removeParameter(), IParameterized::XComponentName(), IParameterized::YComponentName(), and IParameterized::ZComponentName().

Referenced by IParticle::registerPosition().

Here is the call graph for this function:

◆ XComponentName()

std::string IParameterized::XComponentName ( const std::string &  base_name)
staticinherited

◆ YComponentName()

std::string IParameterized::YComponentName ( const std::string &  base_name)
staticinherited

Definition at line 105 of file IParameterized.cpp.

106 {
107  return base_name + "Y";
108 }

Referenced by IParameterized::registerVector(), IParameterized::removeVector(), and IParameterized::setVectorValue().

◆ ZComponentName()

std::string IParameterized::ZComponentName ( const std::string &  base_name)
staticinherited

Definition at line 110 of file IParameterized.cpp.

111 {
112  return base_name + "Z";
113 }

Referenced by IParameterized::registerVector(), IParameterized::removeVector(), and IParameterized::setVectorValue().

◆ setName()

void IParameterized::setName ( const std::string &  name)
inlineinherited

Definition at line 68 of file IParameterized.h.

68 { m_name = name; }
std::string m_name

References IParameterized::m_name.

Referenced by BasicLattice::BasicLattice(), Beam::Beam(), Layer::clone(), ConvolutionDetectorResolution::ConvolutionDetectorResolution(), LayersWithAbsorptionBuilder::createSampleByIndex(), Basic2DParaCrystalBuilder::createSampleByIndex(), ParticleInVacuumBuilder::createSampleByIndex(), SimpleMagneticRotationBuilder::createSampleByIndex(), Crystal::Crystal(), DetectionProperties::DetectionProperties(), DistributionHandler::DistributionHandler(), FormFactorBAPol::FormFactorBAPol(), FormFactorCoreShell::FormFactorCoreShell(), FormFactorCrystal::FormFactorCrystal(), FormFactorDecoratorMaterial::FormFactorDecoratorMaterial(), FormFactorDecoratorPositionFactor::FormFactorDecoratorPositionFactor(), FormFactorDecoratorRotation::FormFactorDecoratorRotation(), FormFactorDWBA::FormFactorDWBA(), FormFactorDWBAPol::FormFactorDWBAPol(), FormFactorWeighted::FormFactorWeighted(), HexagonalLattice::HexagonalLattice(), IDetector::IDetector(), DepthProbeSimulation::initialize(), GISASSimulation::initialize(), OffSpecSimulation::initialize(), SpecularSimulation::initialize(), SpecularDetector1D::initialize(), MesoCrystal::initialize(), Particle::initialize(), ParticleComposition::initialize(), INode::INode(), Instrument::Instrument(), InterferenceFunction1DLattice::InterferenceFunction1DLattice(), InterferenceFunction2DLattice::InterferenceFunction2DLattice(), InterferenceFunction2DParaCrystal::InterferenceFunction2DParaCrystal(), InterferenceFunction2DSuperLattice::InterferenceFunction2DSuperLattice(), InterferenceFunction3DLattice::InterferenceFunction3DLattice(), InterferenceFunctionFinite2DLattice::InterferenceFunctionFinite2DLattice(), InterferenceFunctionFinite3DLattice::InterferenceFunctionFinite3DLattice(), InterferenceFunctionHardDisk::InterferenceFunctionHardDisk(), InterferenceFunctionNone::InterferenceFunctionNone(), InterferenceFunctionRadialParaCrystal::InterferenceFunctionRadialParaCrystal(), InterferenceFunctionTwin::InterferenceFunctionTwin(), ISampleBuilder::ISampleBuilder(), IsGISAXSDetector::IsGISAXSDetector(), Lattice::Lattice(), Layer::Layer(), LayerInterface::LayerInterface(), LayerRoughness::LayerRoughness(), MultiLayer::MultiLayer(), Beam::operator=(), SampleBuilderNode::operator=(), ParticleCoreShell::ParticleCoreShell(), ParticleDistribution::ParticleDistribution(), ParticleLayout::ParticleLayout(), RectangularDetector(), SampleBuilderNode::reset(), ResolutionFunction2DGaussian::ResolutionFunction2DGaussian(), SampleBuilderNode::SampleBuilderNode(), SampleBuilderNode::setSBN(), SphericalDetector::SphericalDetector(), and SquareLattice::SquareLattice().

◆ getName()

Member Data Documentation

◆ m_normal_to_detector

kvector_t RectangularDetector::m_normal_to_detector
private

◆ m_u0

double RectangularDetector::m_u0
private

◆ m_v0

double RectangularDetector::m_v0
private

position of normal vector hitting point in detector coordinates

Definition at line 102 of file RectangularDetector.h.

Referenced by createPixel(), getV0(), indexOfSpecular(), initUandV(), regionOfInterestPixel(), setDistanceAndOffset(), and setPosition().

◆ m_direction

kvector_t RectangularDetector::m_direction
private

direction vector of detector coordinate system

Definition at line 103 of file RectangularDetector.h.

Referenced by getDirectionVector(), initUandV(), and setPosition().

◆ m_distance

double RectangularDetector::m_distance
private

distance from sample origin to the detector plane

Definition at line 104 of file RectangularDetector.h.

Referenced by getDistance(), indexOfSpecular(), initNormalVector(), initUandV(), setDistanceAndOffset(), and setPosition().

◆ m_dbeam_u0

double RectangularDetector::m_dbeam_u0
private

Definition at line 105 of file RectangularDetector.h.

Referenced by getDirectBeamU0(), initUandV(), and setDirectBeamPosition().

◆ m_dbeam_v0

double RectangularDetector::m_dbeam_v0
private

position of direct beam in detector coordinates

Definition at line 105 of file RectangularDetector.h.

Referenced by getDirectBeamV0(), initUandV(), and setDirectBeamPosition().

◆ m_detector_arrangement

◆ m_u_unit

kvector_t RectangularDetector::m_u_unit
private

◆ m_v_unit

kvector_t RectangularDetector::m_v_unit
private

◆ m_detector_mask

◆ m_region_of_interest

std::unique_ptr<RegionOfInterest> IDetector2D::m_region_of_interest
privateinherited

◆ m_axes

◆ m_detection_properties

DetectionProperties IDetector::m_detection_properties
privateinherited

◆ mP_detector_resolution

◆ m_parent

const INode* INode::m_parent {nullptr}
privateinherited

Definition at line 81 of file INode.h.

Referenced by INode::displayName(), INode::parent(), and INode::setParent().

◆ m_NP

const size_t INode::m_NP
protectedinherited

Definition at line 86 of file INode.h.

Referenced by INode::INode().

◆ m_P

std::vector<double> INode::m_P
protectedinherited

Definition at line 87 of file INode.h.

Referenced by INode::INode(), and IFootprintFactor::setWidthRatio().

◆ m_name

std::string IParameterized::m_name
privateinherited

Definition at line 72 of file IParameterized.h.

Referenced by IParameterized::getName(), and IParameterized::setName().

◆ m_pool

std::unique_ptr<ParameterPool> IParameterized::m_pool
privateinherited

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