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

Public Member Functions

 AngularSpecScan (double wl, std::vector< double > inc_angle)
 
 AngularSpecScan (double wl, const IAxis &inc_angle)
 
 AngularSpecScan (double wl, int nbins, double alpha_i_min, double alpha_i_max)
 
 ~AngularSpecScan () override
 
AngularSpecScanclone () const override
 
std::vector< SpecularSimulationElementgenerateSimulationElements () const override
 
virtual const IAxiscoordinateAxis () const override
 
virtual const IFootprintFactorfootprintFactor () const override
 
std::vector< double > footprint (size_t i, size_t n_elements) const override
 
size_t numberOfSimulationElements () const override
 
std::vector< double > createIntensities (const std::vector< SpecularSimulationElement > &sim_elements) const override
 
std::string print () const override
 
double wavelength () const
 
const ScanResolutionwavelengthResolution () const
 
const ScanResolutionangleResolution () const
 
void setFootprintFactor (const IFootprintFactor *f_factor)
 
void setWavelengthResolution (const ScanResolution &resolution)
 
void setRelativeWavelengthResolution (const RangedDistribution &distr, double rel_dev)
 
void setRelativeWavelengthResolution (const RangedDistribution &distr, const std::vector< double > &rel_dev)
 
void setAbsoluteWavelengthResolution (const RangedDistribution &distr, double std_dev)
 
void setAbsoluteWavelengthResolution (const RangedDistribution &distr, const std::vector< double > &std_dev)
 
void setAngleResolution (const ScanResolution &resolution)
 
void setRelativeAngularResolution (const RangedDistribution &distr, double rel_dev)
 
void setRelativeAngularResolution (const RangedDistribution &distr, const std::vector< double > &rel_dev)
 
void setAbsoluteAngularResolution (const RangedDistribution &distr, double std_dev)
 
void setAbsoluteAngularResolution (const RangedDistribution &distr, const std::vector< double > &std_dev)
 
virtual void transferToCPP ()
 

Private Types

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

Private Member Functions

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

Private Attributes

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

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 112 of file AngularSpecScan.h.

Constructor & Destructor Documentation

◆ AngularSpecScan() [1/3]

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

Definition at line 46 of file AngularSpecScan.cpp.

47  : m_wl(wl),
48  m_inc_angle(std::make_unique<PointwiseAxis>("inc_angles", std::move(inc_angle))),
51 {
53 }
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 55 of file AngularSpecScan.cpp.

56  : m_wl(wl), m_inc_angle(inc_angle.clone()),
59 {
61 }
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 63 of file AngularSpecScan.cpp.

64  : m_wl(wl),
65  m_inc_angle(std::make_unique<FixedBinAxis>("inc_angles", nbins, alpha_i_min, alpha_i_max)),
68 {
70 }

References checkInitialization().

Here is the call graph for this function:

◆ ~AngularSpecScan()

AngularSpecScan::~AngularSpecScan ( )
overridedefault

Member Function Documentation

◆ clone()

AngularSpecScan * AngularSpecScan::clone ( ) const
overridevirtual

Implements ISpecularScan.

Definition at line 72 of file AngularSpecScan.cpp.

73 {
74  auto* result = new AngularSpecScan(m_wl, *m_inc_angle);
75  result->setFootprintFactor(m_footprint.get());
76  result->setWavelengthResolution(*m_wl_resolution);
77  result->setAngleResolution(*m_inc_resolution);
78  return result;
79 }
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:

◆ generateSimulationElements()

std::vector< SpecularSimulationElement > AngularSpecScan::generateSimulationElements ( ) 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(
98  SpecularSimulationElement(wl, -inc, wl >= 0 && inc >= 0 && inc <= M_PI_2));
99  }
100  }
101  }
102  return result;
103 }
#define M_PI_2
Definition: MathConstants.h:40
size_t numberOfSimulationElements() const override
Returns the number of simulation elements.
DistrOutput applyIncResolution() const
DistrOutput applyWlResolution() const
A parameter value with a weight, as obtained when sampling from a distribution.
Data stucture containing both input and output of a single image pixel for specular simulation.
std::vector< std::vector< double > > extractValues(std::vector< std::vector< ParameterSample >> samples, const std::function< double(const ParameterSample &)> extractor)

References applyIncResolution(), applyWlResolution(), anonymous_namespace{AngularSpecScan.cpp}::extractValues(), m_inc_angle, M_PI_2, numberOfSimulationElements(), and ParameterSample::value.

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.

Referenced by anonymous_namespace{SpecularSimulation.cpp}::mangledScan(), and print().

◆ 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.

Referenced by anonymous_namespace{SpecularSimulation.cpp}::mangledScan().

◆ 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 }
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...

References applyIncResolution(), anonymous_namespace{AngularSpecScan.cpp}::extractValues(), m_footprint, m_inc_resolution, M_PI_2, m_wl_resolution, 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().

◆ 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].getIntensity() * inc_weight * wl_weights[i][j];
246  ++elem_pos;
247  }
248  }
249  }
250  return result;
251 }

References applyIncResolution(), applyWlResolution(), anonymous_namespace{AngularSpecScan.cpp}::extractValues(), m_inc_angle, and ParameterSample::weight.

Here is the call graph for this function:

◆ print()

std::string AngularSpecScan::print ( ) const
overridevirtual

Print scan definition in python format.

Implements ISpecularScan.

Definition at line 253 of file AngularSpecScan.cpp.

254 {
255  std::stringstream result;
256  result << "\n" << pyfmt::indent() << "# Defining specular scan:\n";
257  const std::string axis_def = pyfmt::indent() + "axis = ";
258  result << axis_def << coordinateAxis()->pyString("rad", axis_def.size()) << "\n";
259 
260  result << pyfmt::indent() << "scan = ";
261  result << "ba.AngularSpecScan(" << pyfmt::printDouble(m_wl) << ", axis)\n";
262 
263  if (m_footprint) {
264  result << *m_footprint << "\n";
265  result << pyfmt::indent() << "scan.setFootprintFactor(footprint)\n";
266  }
267  if (!m_inc_resolution->empty()) {
268  result << "\n" << pyfmt::indent() << "# Defining angular resolution\n";
269  result << *m_inc_resolution << "\n";
270  result << pyfmt::indent() << "scan.setAngleResolution(resolution)\n";
271  }
272  if (!m_wl_resolution->empty()) {
273  result << "\n" << pyfmt::indent() << "# Defining wavelength resolution\n";
274  result << *m_wl_resolution << "\n";
275  result << pyfmt::indent() << "scan.setWavelengthResolution(resolution)\n";
276  }
277  return result.str();
278 }
virtual const IAxis * coordinateAxis() const override
Returns coordinate axis assigned to the data holder.
virtual std::string pyString(const std::string &units, size_t offset) const =0
std::string printDouble(double input)
Definition: PyFmt.cpp:43
std::string indent(size_t width)
Returns a string of blanks with given width.
Definition: PyFmt.cpp:141

References coordinateAxis(), pyfmt::indent(), m_footprint, m_inc_resolution, m_wl, m_wl_resolution, pyfmt::printDouble(), and IAxis::pyString().

Here is the call graph for this function:

◆ wavelength()

double AngularSpecScan::wavelength ( ) const
inline

Definition at line 65 of file AngularSpecScan.h.

65 { return m_wl; }

References m_wl.

◆ wavelengthResolution()

const ScanResolution* AngularSpecScan::wavelengthResolution ( ) const
inline

Definition at line 68 of file AngularSpecScan.h.

68 { return m_wl_resolution.get(); }

References m_wl_resolution.

Referenced by anonymous_namespace{SpecularSimulation.cpp}::mangledScan().

◆ angleResolution()

const ScanResolution* AngularSpecScan::angleResolution ( ) const
inline

Definition at line 69 of file AngularSpecScan.h.

69 { return m_inc_resolution.get(); }

References m_inc_resolution.

Referenced by anonymous_namespace{SpecularSimulation.cpp}::mangledScan().

◆ 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:

◆ 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 }
DistrOutput m_wl_res_cache
ScanResolution * clone() const override=0

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

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

Here is the call graph for this function:

◆ setRelativeWavelengthResolution() [1/2]

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

Definition at line 117 of file AngularSpecScan.cpp.

119 {
120  std::unique_ptr<ScanResolution> resolution(
122  setWavelengthResolution(*resolution);
123 }
void setWavelengthResolution(const ScanResolution &resolution)
Sets wavelength resolution values via ScanResolution object.
static ScanResolution * scanRelativeResolution(const RangedDistribution &distr, double stddev)

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

Here is the call graph for this function:

◆ setRelativeWavelengthResolution() [2/2]

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

Sets wavelength resolution values via RangedDistribution 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:

◆ setAbsoluteWavelengthResolution() [1/2]

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

Definition at line 133 of file AngularSpecScan.cpp.

135 {
136  std::unique_ptr<ScanResolution> resolution(
138  setWavelengthResolution(*resolution);
139 }
static ScanResolution * scanAbsoluteResolution(const RangedDistribution &distr, double stddev)

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

Here is the call graph for this function:

◆ setAbsoluteWavelengthResolution() [2/2]

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

Sets wavelength resolution values via RangedDistribution 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 }

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 }
DistrOutput m_inc_res_cache

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

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

Here is the call graph for this function:

◆ setRelativeAngularResolution() [1/2]

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

Definition at line 156 of file AngularSpecScan.cpp.

157 {
158  std::unique_ptr<ScanResolution> resolution(
160  setAngleResolution(*resolution);
161 }
void setAngleResolution(const ScanResolution &resolution)
Sets angle resolution values via ScanResolution object.

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

Here is the call graph for this function:

◆ setRelativeAngularResolution() [2/2]

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

Sets angular resolution values via RangedDistribution 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 }

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

Here is the call graph for this function:

◆ setAbsoluteAngularResolution() [1/2]

void AngularSpecScan::setAbsoluteAngularResolution ( const RangedDistribution 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:

◆ setAbsoluteAngularResolution() [2/2]

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

Sets angular resolution values via RangedDistribution 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 }

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

Here is the call graph for this function:

◆ checkInitialization()

void AngularSpecScan::checkInitialization ( )
private

Definition at line 280 of file AngularSpecScan.cpp.

281 {
282  if (m_wl <= 0.0)
283  throw std::runtime_error(
284  "Error in AngularSpecScan::checkInitialization: wavelength shell be positive");
285 
286  const std::vector<double> axis_values = m_inc_angle->getBinCenters();
287  if (!std::is_sorted(axis_values.begin(), axis_values.end()))
288  throw std::runtime_error("Error in AngularSpecScan::checkInitialization: q-vector values "
289  "shall be sorted in ascending order.");
290 
291  // TODO: check for inclination angle limits after switching to pointwise resolution.
292 }

References m_inc_angle, and m_wl.

Referenced by AngularSpecScan().

◆ applyWlResolution()

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

Definition at line 294 of file AngularSpecScan.cpp.

295 {
296  if (m_wl_res_cache.empty())
297  m_wl_res_cache = m_wl_resolution->generateSamples(m_wl, m_inc_angle->size());
298  return m_wl_res_cache;
299 }

References m_inc_angle, m_wl, m_wl_res_cache, and m_wl_resolution.

Referenced by createIntensities(), and generateSimulationElements().

◆ applyIncResolution()

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

Definition at line 301 of file AngularSpecScan.cpp.

302 {
303  if (m_inc_res_cache.empty())
304  m_inc_res_cache = m_inc_resolution->generateSamples(m_inc_angle->getBinCenters());
305  return m_inc_res_cache;
306 }

References m_inc_angle, m_inc_res_cache, and m_inc_resolution.

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

◆ transferToCPP()

virtual void ICloneable::transferToCPP ( )
inlinevirtualinherited

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

Definition at line 34 of file ICloneable.h.

Member Data Documentation

◆ m_wl

double AngularSpecScan::m_wl
private

Definition at line 118 of file AngularSpecScan.h.

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

◆ m_inc_angle

◆ m_footprint

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

Definition at line 120 of file AngularSpecScan.h.

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

◆ m_wl_resolution

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

◆ m_wl_res_cache

DistrOutput AngularSpecScan::m_wl_res_cache
mutableprivate

Definition at line 123 of file AngularSpecScan.h.

Referenced by applyWlResolution(), and setWavelengthResolution().

◆ m_inc_resolution

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

◆ m_inc_res_cache

DistrOutput AngularSpecScan::m_inc_res_cache
mutableprivate

Definition at line 126 of file AngularSpecScan.h.

Referenced by applyIncResolution(), and setAngleResolution().


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