BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
AngularSpecScan Class Reference

Scan type with inclination angles as coordinate values and a unique wavelength. More...

Inheritance diagram for AngularSpecScan:
[legend]
Collaboration diagram for AngularSpecScan:
[legend]

Public Member Functions

 AngularSpecScan (double wl, const IAxis &inc_angle)
 
 AngularSpecScan (double wl, int nbins, double alpha_i_min, double alpha_i_max)
 Sets angle-defined specular scan. More...
 
 AngularSpecScan (double wl, std::vector< double > inc_angle)
 
 ~AngularSpecScan () override
 
const ScanResolutionangleResolution () const
 
AngularSpecScanclone () const override
 
virtual const IAxiscoordinateAxis () const override
 Returns coordinate axis assigned to the data holder. More...
 
std::vector< double > createIntensities (const std::vector< SpecularSimulationElement > &sim_elements) const override
 Returns intensity vector corresponding to convolution of given simulation elements. More...
 
std::vector< double > footprint (size_t i, size_t n_elements) const override
 Returns footprint correction factor for a range of simulation elements of size n_elements and starting from element with index i. More...
 
virtual const IFootprintFactorfootprintFactor () const override
 Returns IFootprintFactor object pointer. More...
 
std::vector< SpecularSimulationElementgenerateSimulationElements (const Instrument &instrument) const override
 Generates simulation elements for specular simulations. More...
 
size_t numberOfSimulationElements () const override
 Returns the number of simulation elements. More...
 
void setAbsoluteAngularResolution (const IRangedDistribution &distr, const std::vector< double > &std_dev)
 Sets angular resolution values via IRangedDistribution and values of standard deviations. More...
 
void setAbsoluteAngularResolution (const IRangedDistribution &distr, double std_dev)
 
void setAbsoluteWavelengthResolution (const IRangedDistribution &distr, const std::vector< double > &std_dev)
 Sets wavelength resolution values via IRangedDistribution and values of standard deviations. More...
 
void setAbsoluteWavelengthResolution (const IRangedDistribution &distr, double std_dev)
 
void setAngleResolution (const ScanResolution &resolution)
 Sets angle resolution values via ScanResolution object. More...
 
void setFootprintFactor (const IFootprintFactor *f_factor)
 Sets footprint correction factor. More...
 
void setRelativeAngularResolution (const IRangedDistribution &distr, const std::vector< double > &rel_dev)
 Sets angular resolution values via IRangedDistribution and values of relative deviations (that is, rel_dev equals standard deviation divided by the mean value). More...
 
void setRelativeAngularResolution (const IRangedDistribution &distr, double rel_dev)
 
void setRelativeWavelengthResolution (const IRangedDistribution &distr, const std::vector< double > &rel_dev)
 Sets wavelength resolution values via IRangedDistribution and values of relative deviations (that is, rel_dev equals standard deviation divided by the mean value). More...
 
void setRelativeWavelengthResolution (const IRangedDistribution &distr, double rel_dev)
 
void setWavelengthResolution (const ScanResolution &resolution)
 Sets wavelength resolution values via ScanResolution object. More...
 
virtual void transferToCPP ()
 Used for Python overriding of clone (see swig/tweaks.py) More...
 
double wavelength () const
 
const ScanResolutionwavelengthResolution () const
 

Private Types

using DistrOutput = std::vector< std::vector< ParameterSample > >
 

Private Member Functions

DistrOutput applyIncResolution () const
 
DistrOutput applyWlResolution () const
 
void checkInitialization ()
 

Private Attributes

std::unique_ptr< IFootprintFactorm_footprint
 
const std::unique_ptr< IAxism_inc_angle
 
DistrOutput m_inc_res_cache
 
std::unique_ptr< ScanResolutionm_inc_resolution
 
const double m_wl
 
DistrOutput m_wl_res_cache
 
std::unique_ptr< ScanResolutionm_wl_resolution
 

Detailed Description

Scan type with inclination angles as coordinate values and a unique wavelength.

Features footprint correction.

Definition at line 27 of file AngularSpecScan.h.

Member Typedef Documentation

◆ DistrOutput

using AngularSpecScan::DistrOutput = std::vector<std::vector<ParameterSample> >
private

Definition at line 109 of file AngularSpecScan.h.

Constructor & Destructor Documentation

◆ AngularSpecScan() [1/3]

AngularSpecScan::AngularSpecScan ( double  wl,
std::vector< double >  inc_angle 
)

Definition at line 44 of file AngularSpecScan.cpp.

45  : m_wl(wl)
46  , m_inc_angle(std::make_unique<PointwiseAxis>("inc_angles", std::move(inc_angle)))
49 {
51 }
const double m_wl
const std::unique_ptr< IAxis > m_inc_angle
std::unique_ptr< ScanResolution > m_inc_resolution
std::unique_ptr< ScanResolution > m_wl_resolution
static ScanResolution * scanEmptyResolution()

References checkInitialization().

Referenced by clone().

Here is the call graph for this function:

◆ AngularSpecScan() [2/3]

AngularSpecScan::AngularSpecScan ( double  wl,
const IAxis inc_angle 
)

Definition at line 53 of file AngularSpecScan.cpp.

54  : m_wl(wl)
55  , m_inc_angle(inc_angle.clone())
58 {
60 }
virtual IAxis * clone() const =0
clone function

References checkInitialization().

Here is the call graph for this function:

◆ AngularSpecScan() [3/3]

AngularSpecScan::AngularSpecScan ( double  wl,
int  nbins,
double  alpha_i_min,
double  alpha_i_max 
)

Sets angle-defined specular scan.

The first parameter is always a wavelength in nm. Second parameter is either a numpy array of incident angles in radians or an IAxis object with angle values. Alternatively an axis can be defined in-place, then the second passed parameter is the number of bins, third - minimum on-axis angle value, fourth - maximum on-axis angle value.

Definition at line 62 of file AngularSpecScan.cpp.

63  : m_wl(wl)
64  , m_inc_angle(std::make_unique<FixedBinAxis>("inc_angles", nbins, alpha_i_min, alpha_i_max))
67 {
69 }

References checkInitialization().

Here is the call graph for this function:

◆ ~AngularSpecScan()

AngularSpecScan::~AngularSpecScan ( )
overridedefault

Member Function Documentation

◆ angleResolution()

const ScanResolution* AngularSpecScan::angleResolution ( ) const
inline

Definition at line 66 of file AngularSpecScan.h.

66 { return m_inc_resolution.get(); }

References m_inc_resolution.

Referenced by TransformFromDomain::setSpecularBeamItem().

◆ applyIncResolution()

AngularSpecScan::DistrOutput AngularSpecScan::applyIncResolution ( ) const
private

Definition at line 274 of file AngularSpecScan.cpp.

275 {
276  if (m_inc_res_cache.empty())
277  m_inc_res_cache = m_inc_resolution->generateSamples(m_inc_angle->binCenters());
278  return m_inc_res_cache;
279 }
DistrOutput m_inc_res_cache

References m_inc_angle, m_inc_res_cache, and m_inc_resolution.

Referenced by createIntensities(), footprint(), and generateSimulationElements().

◆ applyWlResolution()

AngularSpecScan::DistrOutput AngularSpecScan::applyWlResolution ( ) const
private

Definition at line 267 of file AngularSpecScan.cpp.

268 {
269  if (m_wl_res_cache.empty())
270  m_wl_res_cache = m_wl_resolution->generateSamples(m_wl, m_inc_angle->size());
271  return m_wl_res_cache;
272 }
DistrOutput m_wl_res_cache

References m_inc_angle, m_wl, m_wl_res_cache, and m_wl_resolution.

Referenced by createIntensities(), and generateSimulationElements().

◆ checkInitialization()

void AngularSpecScan::checkInitialization ( )
private

Definition at line 253 of file AngularSpecScan.cpp.

254 {
255  if (m_wl <= 0.0)
256  throw std::runtime_error(
257  "Error in AngularSpecScan::checkInitialization: wavelength shell be positive");
258 
259  const std::vector<double> axis_values = m_inc_angle->binCenters();
260  if (!std::is_sorted(axis_values.begin(), axis_values.end()))
261  throw std::runtime_error("Error in AngularSpecScan::checkInitialization: q-vector values "
262  "shall be sorted in ascending order.");
263 
264  // TODO: check for inclination angle limits after switching to pointwise resolution.
265 }

References m_inc_angle, and m_wl.

Referenced by AngularSpecScan().

◆ clone()

AngularSpecScan * AngularSpecScan::clone ( ) const
overridevirtual

Implements ISpecularScan.

Definition at line 71 of file AngularSpecScan.cpp.

72 {
73  auto* result = new AngularSpecScan(m_wl, *m_inc_angle);
74  result->setFootprintFactor(m_footprint.get());
75  result->setWavelengthResolution(*m_wl_resolution);
76  result->setAngleResolution(*m_inc_resolution);
77  return result;
78 }
AngularSpecScan(double wl, std::vector< double > inc_angle)
std::unique_ptr< IFootprintFactor > m_footprint

References AngularSpecScan(), m_footprint, m_inc_angle, m_inc_resolution, m_wl, and m_wl_resolution.

Here is the call graph for this function:

◆ coordinateAxis()

virtual const IAxis* AngularSpecScan::coordinateAxis ( ) const
inlineoverridevirtual

Returns coordinate axis assigned to the data holder.

Implements ISpecularScan.

Definition at line 46 of file AngularSpecScan.h.

46 { return m_inc_angle.get(); }

References m_inc_angle.

◆ createIntensities()

std::vector< double > AngularSpecScan::createIntensities ( const std::vector< SpecularSimulationElement > &  sim_elements) const
overridevirtual

Returns intensity vector corresponding to convolution of given simulation elements.

Implements ISpecularScan.

Definition at line 229 of file AngularSpecScan.cpp.

230 {
231  const size_t axis_size = m_inc_angle->size();
232  std::vector<double> result(axis_size, 0.0);
233 
234  const auto wl_weights = extractValues(
235  applyWlResolution(), [](const ParameterSample& sample) { return sample.weight; });
236  const auto inc_weights = extractValues(
237  applyIncResolution(), [](const ParameterSample& sample) { return sample.weight; });
238 
239  size_t elem_pos = 0;
240  for (size_t i = 0; i < axis_size; ++i) {
241  double& current = result[i];
242  for (size_t k = 0, size_incs = inc_weights[i].size(); k < size_incs; ++k) {
243  const double inc_weight = inc_weights[i][k];
244  for (size_t j = 0, size_wls = wl_weights[i].size(); j < size_wls; ++j) {
245  current += sim_elements[elem_pos].intensity() * inc_weight * wl_weights[i][j];
246  ++elem_pos;
247  }
248  }
249  }
250  return result;
251 }
DistrOutput applyIncResolution() const
DistrOutput applyWlResolution() const
A parameter value with a weight, as obtained when sampling from a distribution.

References applyIncResolution(), applyWlResolution(), m_inc_angle, and ParameterSample::weight.

Here is the call graph for this function:

◆ footprint()

std::vector< double > AngularSpecScan::footprint ( size_t  i,
size_t  n_elements 
) const
overridevirtual

Returns footprint correction factor for a range of simulation elements of size n_elements and starting from element with index i.

Implements ISpecularScan.

Definition at line 186 of file AngularSpecScan.cpp.

187 {
188  if (start + n_elements > numberOfSimulationElements())
189  throw std::runtime_error("Error in AngularSpecScan::footprint: given index exceeds the "
190  "number of simulation elements");
191 
192  std::vector<double> result(n_elements, 1.0);
193  if (!m_footprint)
194  return result;
195 
196  const size_t n_wl_samples = m_wl_resolution->nSamples();
197  const size_t n_inc_samples = m_inc_resolution->nSamples();
198 
199  const auto sample_values = extractValues(
200  applyIncResolution(), [](const ParameterSample& sample) { return sample.value; });
201 
202  const size_t pos_out = start / (n_wl_samples * n_inc_samples);
203  size_t pos_inc = (start - pos_out * n_wl_samples * n_inc_samples) / n_wl_samples;
204  size_t pos_wl = (start - pos_inc * n_wl_samples);
205  int left = static_cast<int>(n_elements);
206  size_t pos_res = 0;
207  for (size_t i = pos_out; left > 0; ++i)
208  for (size_t k = pos_inc; k < n_inc_samples && left > 0; ++k) {
209  pos_inc = 0;
210  const double angle = sample_values[i][k];
211  const double footprint =
212  (angle >= 0 && angle <= M_PI_2) ? m_footprint->calculate(angle) : 1.0;
213  for (size_t j = pos_wl; j < n_wl_samples && left > 0; ++j) {
214  pos_wl = 0;
215  result[pos_res] = footprint;
216  ++pos_res;
217  --left;
218  }
219  }
220  return result;
221 }
#define M_PI_2
Definition: Constants.h:45
std::vector< double > footprint(size_t i, size_t n_elements) const override
Returns footprint correction factor for a range of simulation elements of size n_elements and startin...
size_t numberOfSimulationElements() const override
Returns the number of simulation elements.

References applyIncResolution(), m_footprint, m_inc_resolution, M_PI_2, m_wl_resolution, numberOfSimulationElements(), and ParameterSample::value.

Here is the call graph for this function:

◆ footprintFactor()

virtual const IFootprintFactor* AngularSpecScan::footprintFactor ( ) const
inlineoverridevirtual

Returns IFootprintFactor object pointer.

Implements ISpecularScan.

Definition at line 49 of file AngularSpecScan.h.

49 { return m_footprint.get(); }

References m_footprint.

◆ generateSimulationElements()

std::vector< SpecularSimulationElement > AngularSpecScan::generateSimulationElements ( const Instrument instrument) const
overridevirtual

Generates simulation elements for specular simulations.

Implements ISpecularScan.

Definition at line 83 of file AngularSpecScan.cpp.

84 {
85  const auto wls = extractValues(applyWlResolution(),
86  [](const ParameterSample& sample) { return sample.value; });
87  const auto incs = extractValues(applyIncResolution(),
88  [](const ParameterSample& sample) { return sample.value; });
89 
90  std::vector<SpecularSimulationElement> result;
91  result.reserve(numberOfSimulationElements());
92  for (size_t i = 0; i < m_inc_angle->size(); ++i) {
93  for (size_t k = 0, size_incs = incs[i].size(); k < size_incs; ++k) {
94  const double inc = incs[i][k];
95  for (size_t j = 0, size_wls = wls[i].size(); j < size_wls; ++j) {
96  const double wl = wls[i][j];
97  result.emplace_back(SpecularSimulationElement(
98  wl, -inc, instrument, wl >= 0 && inc >= 0 && inc <= M_PI_2));
99  }
100  }
101  }
102  return result;
103 }
Data stucture containing both input and output of a single image pixel for specular simulation.

References applyIncResolution(), applyWlResolution(), m_inc_angle, M_PI_2, numberOfSimulationElements(), and ParameterSample::value.

Here is the call graph for this function:

◆ numberOfSimulationElements()

size_t AngularSpecScan::numberOfSimulationElements ( ) const
overridevirtual

Returns the number of simulation elements.

Implements ISpecularScan.

Definition at line 223 of file AngularSpecScan.cpp.

224 {
225  return m_inc_angle->size() * m_wl_resolution->nSamples() * m_inc_resolution->nSamples();
226 }

References m_inc_angle, m_inc_resolution, and m_wl_resolution.

Referenced by footprint(), and generateSimulationElements().

◆ setAbsoluteAngularResolution() [1/2]

void AngularSpecScan::setAbsoluteAngularResolution ( const IRangedDistribution distr,
const std::vector< double > &  std_dev 
)

Sets angular resolution values via IRangedDistribution and values of standard deviations.

std_dev can be either single-valued or a numpy array. In the latter case the length of the array should coinside with the length of the inclination angle axis.

Definition at line 178 of file AngularSpecScan.cpp.

180 {
181  std::unique_ptr<ScanResolution> resolution(
183  setAngleResolution(*resolution);
184 }
void setAngleResolution(const ScanResolution &resolution)
Sets angle resolution values via ScanResolution object.
static ScanResolution * scanAbsoluteResolution(const IRangedDistribution &distr, double stddev)

References ScanResolution::scanAbsoluteResolution(), and setAngleResolution().

Here is the call graph for this function:

◆ setAbsoluteAngularResolution() [2/2]

void AngularSpecScan::setAbsoluteAngularResolution ( const IRangedDistribution distr,
double  std_dev 
)

Definition at line 171 of file AngularSpecScan.cpp.

172 {
173  std::unique_ptr<ScanResolution> resolution(
175  setAngleResolution(*resolution);
176 }

References ScanResolution::scanAbsoluteResolution(), and setAngleResolution().

Here is the call graph for this function:

◆ setAbsoluteWavelengthResolution() [1/2]

void AngularSpecScan::setAbsoluteWavelengthResolution ( const IRangedDistribution distr,
const std::vector< double > &  std_dev 
)

Sets wavelength resolution values via IRangedDistribution and values of standard deviations.

std_dev can be either single-valued or a numpy array. In the latter case the length of the array should coinside with the length of the inclination angle axis.

Definition at line 141 of file AngularSpecScan.cpp.

143 {
144  std::unique_ptr<ScanResolution> resolution(
146  setWavelengthResolution(*resolution);
147 }
void setWavelengthResolution(const ScanResolution &resolution)
Sets wavelength resolution values via ScanResolution object.

References ScanResolution::scanAbsoluteResolution(), and setWavelengthResolution().

Here is the call graph for this function:

◆ setAbsoluteWavelengthResolution() [2/2]

void AngularSpecScan::setAbsoluteWavelengthResolution ( const IRangedDistribution distr,
double  std_dev 
)

Definition at line 133 of file AngularSpecScan.cpp.

135 {
136  std::unique_ptr<ScanResolution> resolution(
138  setWavelengthResolution(*resolution);
139 }

References ScanResolution::scanAbsoluteResolution(), and setWavelengthResolution().

Here is the call graph for this function:

◆ setAngleResolution()

void AngularSpecScan::setAngleResolution ( const ScanResolution resolution)

Sets angle resolution values via ScanResolution object.

Definition at line 149 of file AngularSpecScan.cpp.

150 {
151  m_inc_resolution.reset(resolution.clone());
152  m_inc_res_cache.clear();
153  m_inc_res_cache.shrink_to_fit();
154 }
ScanResolution * clone() const override=0

References ScanResolution::clone(), m_inc_res_cache, and m_inc_resolution.

Referenced by TransformToDomain::addBeamDivergencesToScan(), setAbsoluteAngularResolution(), setRelativeAngularResolution(), and StandardSimulations::SpecularDivergentBeam().

Here is the call graph for this function:

◆ setFootprintFactor()

void AngularSpecScan::setFootprintFactor ( const IFootprintFactor f_factor)

Sets footprint correction factor.

Definition at line 105 of file AngularSpecScan.cpp.

106 {
107  m_footprint.reset(f_factor ? f_factor->clone() : nullptr);
108 }
virtual IFootprintFactor * clone() const =0

References IFootprintFactor::clone(), and m_footprint.

Referenced by StandardSimulations::SpecularWithGaussianBeam(), and StandardSimulations::SpecularWithSquareBeam().

Here is the call graph for this function:

◆ setRelativeAngularResolution() [1/2]

void AngularSpecScan::setRelativeAngularResolution ( const IRangedDistribution distr,
const std::vector< double > &  rel_dev 
)

Sets angular resolution values via IRangedDistribution and values of relative deviations (that is, rel_dev equals standard deviation divided by the mean value).

rel_dev can be either single-valued or a numpy array. In the latter case the length of the array should coinside with the length of the inclination angle axis.

Definition at line 163 of file AngularSpecScan.cpp.

165 {
166  std::unique_ptr<ScanResolution> resolution(
168  setAngleResolution(*resolution);
169 }
static ScanResolution * scanRelativeResolution(const IRangedDistribution &distr, double stddev)

References ScanResolution::scanRelativeResolution(), and setAngleResolution().

Here is the call graph for this function:

◆ setRelativeAngularResolution() [2/2]

void AngularSpecScan::setRelativeAngularResolution ( const IRangedDistribution distr,
double  rel_dev 
)

Definition at line 156 of file AngularSpecScan.cpp.

157 {
158  std::unique_ptr<ScanResolution> resolution(
160  setAngleResolution(*resolution);
161 }

References ScanResolution::scanRelativeResolution(), and setAngleResolution().

Here is the call graph for this function:

◆ setRelativeWavelengthResolution() [1/2]

void AngularSpecScan::setRelativeWavelengthResolution ( const IRangedDistribution distr,
const std::vector< double > &  rel_dev 
)

Sets wavelength resolution values via IRangedDistribution and values of relative deviations (that is, rel_dev equals standard deviation divided by the mean value).

rel_dev can be either single-valued or a numpy array. In the latter case the length of the array should coinside with the length of the inclination angle axis.

Definition at line 125 of file AngularSpecScan.cpp.

127 {
128  std::unique_ptr<ScanResolution> resolution(
130  setWavelengthResolution(*resolution);
131 }

References ScanResolution::scanRelativeResolution(), and setWavelengthResolution().

Here is the call graph for this function:

◆ setRelativeWavelengthResolution() [2/2]

void AngularSpecScan::setRelativeWavelengthResolution ( const IRangedDistribution distr,
double  rel_dev 
)

Definition at line 117 of file AngularSpecScan.cpp.

119 {
120  std::unique_ptr<ScanResolution> resolution(
122  setWavelengthResolution(*resolution);
123 }

References ScanResolution::scanRelativeResolution(), and setWavelengthResolution().

Here is the call graph for this function:

◆ setWavelengthResolution()

void AngularSpecScan::setWavelengthResolution ( const ScanResolution resolution)

Sets wavelength resolution values via ScanResolution object.

Definition at line 110 of file AngularSpecScan.cpp.

111 {
112  m_wl_resolution.reset(resolution.clone());
113  m_wl_res_cache.clear();
114  m_wl_res_cache.shrink_to_fit();
115 }

References ScanResolution::clone(), m_wl_res_cache, and m_wl_resolution.

Referenced by TransformToDomain::addBeamDivergencesToScan(), setAbsoluteWavelengthResolution(), setRelativeWavelengthResolution(), and StandardSimulations::SpecularDivergentBeam().

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.

◆ wavelength()

double AngularSpecScan::wavelength ( ) const
inline

Definition at line 62 of file AngularSpecScan.h.

62 { return m_wl; }

References m_wl.

Referenced by TransformFromDomain::setSpecularBeamItem().

◆ wavelengthResolution()

const ScanResolution* AngularSpecScan::wavelengthResolution ( ) const
inline

Definition at line 65 of file AngularSpecScan.h.

65 { return m_wl_resolution.get(); }

References m_wl_resolution.

Referenced by TransformFromDomain::setSpecularBeamItem().

Member Data Documentation

◆ m_footprint

std::unique_ptr<IFootprintFactor> AngularSpecScan::m_footprint
private

Definition at line 117 of file AngularSpecScan.h.

Referenced by clone(), footprint(), footprintFactor(), and setFootprintFactor().

◆ m_inc_angle

const std::unique_ptr<IAxis> AngularSpecScan::m_inc_angle
private

◆ m_inc_res_cache

DistrOutput AngularSpecScan::m_inc_res_cache
mutableprivate

Definition at line 123 of file AngularSpecScan.h.

Referenced by applyIncResolution(), and setAngleResolution().

◆ m_inc_resolution

std::unique_ptr<ScanResolution> AngularSpecScan::m_inc_resolution
private

◆ m_wl

const double AngularSpecScan::m_wl
private

Definition at line 115 of file AngularSpecScan.h.

Referenced by applyWlResolution(), checkInitialization(), clone(), and wavelength().

◆ m_wl_res_cache

DistrOutput AngularSpecScan::m_wl_res_cache
mutableprivate

Definition at line 120 of file AngularSpecScan.h.

Referenced by applyWlResolution(), and setWavelengthResolution().

◆ m_wl_resolution

std::unique_ptr<ScanResolution> AngularSpecScan::m_wl_resolution
private

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