BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
AlphaScan Class Reference

Description

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

Definition at line 27 of file AlphaScan.h.

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

Public Member Functions

 AlphaScan (double wl, const IAxis &inc_angle)
 
 AlphaScan (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. More...
 
 AlphaScan (double wl, std::vector< double > inc_angle)
 
 ~AlphaScan () override
 
const PolFilteranalyzer () const
 
const ScanResolutionangleResolution () const
 
AlphaScanclone () const override
 
const IAxiscoordinateAxis () const override
 Returns coordinate axis assigned to the data holder. More...
 
CoordSystem1DcreateCoordSystem () const override
 
std::vector< double > createIntensities (const std::vector< SpecularElement > &eles) const override
 Returns intensity vector corresponding to convolution of given simulation elements. More...
 
std::vector< double > footprint (size_t start, 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...
 
const IFootprintFactorfootprintFactor () const override
 Returns IFootprintFactor object pointer. More...
 
std::vector< SpecularElementgenerateElements () const override
 Generates simulation elements for specular simulations. More...
 
size_t numberOfElements () const override
 Returns the number of simulation elements. More...
 
bool polarized () const
 
PolMatrices polMatrices () const
 
void 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. 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. 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. More...
 
void setAbsoluteWavelengthResolution (const IRangedDistribution &distr, double std_dev)
 
void setAnalyzer (R3 direction, double efficiency, double total_transmission)
 Sets the polarization analyzer characteristics of the detector. More...
 
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 setPolarization (R3 bloch_vector)
 Sets the polarization density matrix according to the given Bloch vector. 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). 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. 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). 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. 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 override
 
const ScanResolutionwavelengthResolution () const
 

Protected Attributes

std::unique_ptr< R3 > m_beamPolarization
 Bloch vector encoding the beam's polarization. More...
 
std::unique_ptr< PolFilterm_polAnalyzer
 

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
 
std::unique_ptr< ScanResolutionm_inc_resolution
 
const double m_wl
 
std::unique_ptr< ScanResolutionm_wl_resolution
 

Member Typedef Documentation

◆ DistrOutput

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

Definition at line 109 of file AlphaScan.h.

Constructor & Destructor Documentation

◆ AlphaScan() [1/3]

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

Definition at line 48 of file AlphaScan.cpp.

49  : m_wl(wl)
50  , m_inc_angle(std::make_unique<PointwiseAxis>("inc_angles", std::move(inc_angle)))
53 {
55 }
std::unique_ptr< ScanResolution > m_wl_resolution
Definition: AlphaScan.h:119
void checkInitialization()
Definition: AlphaScan.cpp:259
std::unique_ptr< ScanResolution > m_inc_resolution
Definition: AlphaScan.h:120
const double m_wl
Definition: AlphaScan.h:115
const std::unique_ptr< IAxis > m_inc_angle
Definition: AlphaScan.h:116
static ScanResolution * scanEmptyResolution()

References checkInitialization().

Referenced by clone().

Here is the call graph for this function:

◆ AlphaScan() [2/3]

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

Definition at line 57 of file AlphaScan.cpp.

58  : m_wl(wl)
59  , m_inc_angle(inc_angle.clone())
62 {
64 }
virtual IAxis * clone() const =0

References checkInitialization().

Here is the call graph for this function:

◆ AlphaScan() [3/3]

AlphaScan::AlphaScan ( 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 66 of file AlphaScan.cpp.

67  : m_wl(wl)
68  , m_inc_angle(std::make_unique<FixedBinAxis>("inc_angles", nbins, alpha_i_min, alpha_i_max))
71 {
73 }

References checkInitialization().

Here is the call graph for this function:

◆ ~AlphaScan()

AlphaScan::~AlphaScan ( )
overridedefault

Member Function Documentation

◆ analyzer()

const PolFilter* ISpecularScan::analyzer ( ) const
inlineinherited

Definition at line 45 of file ISpecularScan.h.

45 { return m_polAnalyzer.get(); }
std::unique_ptr< PolFilter > m_polAnalyzer
Definition: ISpecularScan.h:78

References ISpecularScan::m_polAnalyzer.

◆ angleResolution()

const ScanResolution* AlphaScan::angleResolution ( ) const
inline

Definition at line 64 of file AlphaScan.h.

64 { return m_inc_resolution.get(); }

References m_inc_resolution.

◆ applyIncResolution()

AlphaScan::DistrOutput AlphaScan::applyIncResolution ( ) const
private

Definition at line 278 of file AlphaScan.cpp.

279 {
280  return m_inc_resolution->generateSamples(m_inc_angle->binCenters());
281 }

References m_inc_angle, and m_inc_resolution.

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

◆ applyWlResolution()

AlphaScan::DistrOutput AlphaScan::applyWlResolution ( ) const
private

Definition at line 273 of file AlphaScan.cpp.

274 {
275  return m_wl_resolution->generateSamples(m_wl, m_inc_angle->size());
276 }

References m_inc_angle, m_wl, and m_wl_resolution.

Referenced by createIntensities(), and generateElements().

◆ checkInitialization()

void AlphaScan::checkInitialization ( )
private

Definition at line 259 of file AlphaScan.cpp.

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

References m_inc_angle, and m_wl.

Referenced by AlphaScan().

◆ clone()

AlphaScan * AlphaScan::clone ( ) const
overridevirtual

Implements ISpecularScan.

Definition at line 75 of file AlphaScan.cpp.

76 {
77  auto* result = new AlphaScan(m_wl, *m_inc_angle);
78  result->setFootprintFactor(m_footprint.get());
79  result->setWavelengthResolution(*m_wl_resolution);
80  result->setAngleResolution(*m_inc_resolution);
82  result->m_beamPolarization.reset(new R3(*m_beamPolarization));
83  if (m_polAnalyzer)
84  result->m_polAnalyzer.reset(new PolFilter(*m_polAnalyzer));
85  return result;
86 }
std::unique_ptr< IFootprintFactor > m_footprint
Definition: AlphaScan.h:117
AlphaScan(double wl, std::vector< double > inc_angle)
Definition: AlphaScan.cpp:48
std::unique_ptr< R3 > m_beamPolarization
Bloch vector encoding the beam's polarization.
Definition: ISpecularScan.h:77
Detector properties (efficiency, transmission).
Definition: PolFilter.h:27

References AlphaScan(), ISpecularScan::m_beamPolarization, m_footprint, m_inc_angle, m_inc_resolution, ISpecularScan::m_polAnalyzer, m_wl, and m_wl_resolution.

Here is the call graph for this function:

◆ coordinateAxis()

const IAxis* AlphaScan::coordinateAxis ( ) const
inlineoverridevirtual

Returns coordinate axis assigned to the data holder.

Implements ISpecularScan.

Definition at line 45 of file AlphaScan.h.

45 { return m_inc_angle.get(); }

References m_inc_angle.

Referenced by createCoordSystem().

◆ createCoordSystem()

CoordSystem1D * AlphaScan::createCoordSystem ( ) const
overridevirtual

Implements ISpecularScan.

Definition at line 230 of file AlphaScan.cpp.

231 {
233 }
const IAxis * coordinateAxis() const override
Returns coordinate axis assigned to the data holder.
Definition: AlphaScan.h:45
double wavelength() const override
Definition: AlphaScan.h:67
Conversion of axis units for the case of conventional (angle-based) reflectometry.
Definition: CoordSystem1D.h:69

References coordinateAxis(), and wavelength().

Here is the call graph for this function:

◆ createIntensities()

std::vector< double > AlphaScan::createIntensities ( const std::vector< SpecularElement > &  eles) const
overridevirtual

Returns intensity vector corresponding to convolution of given simulation elements.

Implements ISpecularScan.

Definition at line 235 of file AlphaScan.cpp.

236 {
237  const size_t axis_size = m_inc_angle->size();
238  std::vector<double> result(axis_size, 0.0);
239 
240  const auto wl_weights = extractValues(
241  applyWlResolution(), [](const ParameterSample& sample) { return sample.weight; });
242  const auto inc_weights = extractValues(
243  applyIncResolution(), [](const ParameterSample& sample) { return sample.weight; });
244 
245  size_t elem_pos = 0;
246  for (size_t i = 0; i < axis_size; ++i) {
247  double& current = result[i];
248  for (size_t k = 0, size_incs = inc_weights[i].size(); k < size_incs; ++k) {
249  const double inc_weight = inc_weights[i][k];
250  for (size_t j = 0, size_wls = wl_weights[i].size(); j < size_wls; ++j) {
251  current += eles[elem_pos].intensity() * inc_weight * wl_weights[i][j];
252  ++elem_pos;
253  }
254  }
255  }
256  return result;
257 }
DistrOutput applyIncResolution() const
Definition: AlphaScan.cpp:278
DistrOutput applyWlResolution() const
Definition: AlphaScan.cpp:273
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 > AlphaScan::footprint ( size_t  start,
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 188 of file AlphaScan.cpp.

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

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

Here is the call graph for this function:

◆ footprintFactor()

const IFootprintFactor* AlphaScan::footprintFactor ( ) const
inlineoverridevirtual

Returns IFootprintFactor object pointer.

Implements ISpecularScan.

Definition at line 48 of file AlphaScan.h.

48 { return m_footprint.get(); }

References m_footprint.

◆ generateElements()

std::vector< SpecularElement > AlphaScan::generateElements ( ) const
overridevirtual

Generates simulation elements for specular simulations.

Implements ISpecularScan.

Definition at line 90 of file AlphaScan.cpp.

91 {
92  const auto wls = extractValues(applyWlResolution(),
93  [](const ParameterSample& sample) { return sample.value; });
94  const auto incs = extractValues(applyIncResolution(),
95  [](const ParameterSample& sample) { return sample.value; });
96 
97  std::vector<SpecularElement> result;
98  result.reserve(numberOfElements());
99  for (size_t i = 0; i < m_inc_angle->size(); ++i) {
100  for (size_t k = 0, size_incs = incs[i].size(); k < size_incs; ++k) {
101  const double alpha = incs[i][k];
102  for (size_t j = 0, size_wls = wls[i].size(); j < size_wls; ++j) {
103  const double wavelength = wls[i][j];
104  const bool computable = wavelength >= 0 && alpha >= 0 && alpha <= M_PI_2;
105  result.emplace_back(
106  SpecularElement::FromAlphaScan(wavelength, -alpha, polMatrices(), computable));
107  }
108  }
109  }
110  return result;
111 }
PolMatrices polMatrices() const
static SpecularElement FromAlphaScan(double wavelength, double alpha, const PolMatrices &polMatrices, bool computable)

References applyIncResolution(), applyWlResolution(), SpecularElement::FromAlphaScan(), m_inc_angle, M_PI_2, numberOfElements(), ISpecularScan::polMatrices(), ParameterSample::value, and wavelength().

Here is the call graph for this function:

◆ numberOfElements()

size_t AlphaScan::numberOfElements ( ) const
overridevirtual

Returns the number of simulation elements.

Implements ISpecularScan.

Definition at line 225 of file AlphaScan.cpp.

226 {
227  return m_inc_angle->size() * m_wl_resolution->nSamples() * m_inc_resolution->nSamples();
228 }

References m_inc_angle, m_inc_resolution, and m_wl_resolution.

Referenced by footprint(), and generateElements().

◆ polarized()

bool ISpecularScan::polarized ( ) const
inherited

Definition at line 31 of file ISpecularScan.cpp.

32 {
34 }

References ISpecularScan::m_beamPolarization, and ISpecularScan::m_polAnalyzer.

◆ polMatrices()

PolMatrices ISpecularScan::polMatrices ( ) const
inherited

Definition at line 36 of file ISpecularScan.cpp.

37 {
38  PolMatrices result;
41  if (m_polAnalyzer)
42  result.setAnalyzerMatrix(m_polAnalyzer->matrix());
43  return result;
44 }
Convenience class for handling polarization density matrix and polarization analyzer operator.
Definition: PolMatrices.h:29
void setAnalyzerMatrix(const SpinMatrix &analyzer)
Sets the polarization analyzer operator (in spin basis along z-axis)
Definition: PolMatrices.h:41
void setPolarizerMatrix(const SpinMatrix &polarization)
Sets the polarization density matrix (in spin basis along z-axis)
Definition: PolMatrices.h:35
static SpinMatrix FromBlochVector(const R3 &v)
Constructs matrix (I+v*s)/2, where s is the Pauli vector.
Definition: SpinMatrix.cpp:41

References SpinMatrix::FromBlochVector(), ISpecularScan::m_beamPolarization, ISpecularScan::m_polAnalyzer, PolMatrices::setAnalyzerMatrix(), and PolMatrices::setPolarizerMatrix().

Referenced by generateElements(), and QzScan::generateElements().

Here is the call graph for this function:

◆ setAbsoluteAngularResolution() [1/2]

void AlphaScan::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 180 of file AlphaScan.cpp.

182 {
183  std::unique_ptr<ScanResolution> resolution(
185  setAngleResolution(*resolution);
186 }
void setAngleResolution(const ScanResolution &resolution)
Sets angle resolution values via ScanResolution object.
Definition: AlphaScan.cpp:153
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 AlphaScan::setAbsoluteAngularResolution ( const IRangedDistribution distr,
double  std_dev 
)

Definition at line 173 of file AlphaScan.cpp.

174 {
175  std::unique_ptr<ScanResolution> resolution(
177  setAngleResolution(*resolution);
178 }

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

Here is the call graph for this function:

◆ setAbsoluteWavelengthResolution() [1/2]

void AlphaScan::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 145 of file AlphaScan.cpp.

147 {
148  std::unique_ptr<ScanResolution> resolution(
150  setWavelengthResolution(*resolution);
151 }
void setWavelengthResolution(const ScanResolution &resolution)
Sets wavelength resolution values via ScanResolution object.
Definition: AlphaScan.cpp:118

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

Here is the call graph for this function:

◆ setAbsoluteWavelengthResolution() [2/2]

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

Definition at line 138 of file AlphaScan.cpp.

139 {
140  std::unique_ptr<ScanResolution> resolution(
142  setWavelengthResolution(*resolution);
143 }

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

Here is the call graph for this function:

◆ setAnalyzer()

void ISpecularScan::setAnalyzer ( R3  direction,
double  efficiency,
double  total_transmission 
)
inherited

Sets the polarization analyzer characteristics of the detector.

Definition at line 26 of file ISpecularScan.cpp.

27 {
28  m_polAnalyzer.reset(new PolFilter(direction, efficiency, total_transmission));
29 }

References ISpecularScan::m_polAnalyzer.

◆ setAngleResolution()

void AlphaScan::setAngleResolution ( const ScanResolution resolution)

Sets angle resolution values via ScanResolution object.

Definition at line 153 of file AlphaScan.cpp.

154 {
155  m_inc_resolution.reset(resolution.clone());
156 }
ScanResolution * clone() const override=0

References ScanResolution::clone(), and m_inc_resolution.

Referenced by setAbsoluteAngularResolution(), and setRelativeAngularResolution().

Here is the call graph for this function:

◆ setFootprintFactor()

void AlphaScan::setFootprintFactor ( const IFootprintFactor f_factor)

Sets footprint correction factor.

Definition at line 113 of file AlphaScan.cpp.

114 {
115  m_footprint.reset(f_factor ? f_factor->clone() : nullptr);
116 }
IFootprintFactor * clone() const override=0

References IFootprintFactor::clone(), and m_footprint.

Here is the call graph for this function:

◆ setPolarization()

void ISpecularScan::setPolarization ( R3  bloch_vector)
inherited

Sets the polarization density matrix according to the given Bloch vector.

Definition at line 21 of file ISpecularScan.cpp.

22 {
23  m_beamPolarization.reset(new R3(bloch_vector));
24 }

References ISpecularScan::m_beamPolarization.

◆ setRelativeAngularResolution() [1/2]

void AlphaScan::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 165 of file AlphaScan.cpp.

167 {
168  std::unique_ptr<ScanResolution> resolution(
170  setAngleResolution(*resolution);
171 }
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 AlphaScan::setRelativeAngularResolution ( const IRangedDistribution distr,
double  rel_dev 
)

Definition at line 158 of file AlphaScan.cpp.

159 {
160  std::unique_ptr<ScanResolution> resolution(
162  setAngleResolution(*resolution);
163 }

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

Here is the call graph for this function:

◆ setRelativeWavelengthResolution() [1/2]

void AlphaScan::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 130 of file AlphaScan.cpp.

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

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

Here is the call graph for this function:

◆ setRelativeWavelengthResolution() [2/2]

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

Definition at line 123 of file AlphaScan.cpp.

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

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

Here is the call graph for this function:

◆ setWavelengthResolution()

void AlphaScan::setWavelengthResolution ( const ScanResolution resolution)

Sets wavelength resolution values via ScanResolution object.

Definition at line 118 of file AlphaScan.cpp.

119 {
120  m_wl_resolution.reset(resolution.clone());
121 }

References ScanResolution::clone(), and m_wl_resolution.

Referenced by setAbsoluteWavelengthResolution(), and setRelativeWavelengthResolution().

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 32 of file ICloneable.h.

◆ wavelength()

double AlphaScan::wavelength ( ) const
inlineoverridevirtual

Implements ISpecularScan.

Definition at line 67 of file AlphaScan.h.

67 { return m_wl; }

References m_wl.

Referenced by createCoordSystem(), and generateElements().

◆ wavelengthResolution()

const ScanResolution* AlphaScan::wavelengthResolution ( ) const
inline

Definition at line 63 of file AlphaScan.h.

63 { return m_wl_resolution.get(); }

References m_wl_resolution.

Member Data Documentation

◆ m_beamPolarization

std::unique_ptr<R3> ISpecularScan::m_beamPolarization
protectedinherited

Bloch vector encoding the beam's polarization.

Definition at line 77 of file ISpecularScan.h.

Referenced by clone(), QzScan::clone(), ISpecularScan::polarized(), ISpecularScan::polMatrices(), and ISpecularScan::setPolarization().

◆ m_footprint

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

Definition at line 117 of file AlphaScan.h.

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

◆ m_inc_angle

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

◆ m_inc_resolution

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

◆ m_polAnalyzer

std::unique_ptr<PolFilter> ISpecularScan::m_polAnalyzer
protectedinherited

◆ m_wl

const double AlphaScan::m_wl
private

Definition at line 115 of file AlphaScan.h.

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

◆ m_wl_resolution

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

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