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

Description

Simulation of radiation depth profile.

Holds an instrument and sample model. Computes radiation intensity as function of incoming glancing angle and penetration depth. Scattered rays are neglected. Only refraction, reflection and attenuation of the incoming beam are accounted for.

Definition at line 37 of file DepthProbeSimulation.h.

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

Public Member Functions

 DepthProbeSimulation (const MultiLayer &sample)
 
 ~DepthProbeSimulation () override
 
void addMask (const IShape2D &shape, bool mask_value=true)
 Adds mask of given shape to the stack of detector masks. The mask value 'true' means that the channel will be excluded from the simulation. The mask which is added last has priority. More...
 
void addParameterDistribution (const ParameterDistribution &par_distr)
 
void addParameterDistribution (ParameterDistribution::WhichParameter whichParameter, const IDistribution1D &distribution, size_t nbr_samples, double sigma_factor=0.0, const RealLimits &limits=RealLimits())
 
const IAxisalphaAxis () const
 Returns a pointer to incident angle axis. More...
 
const IBackgroundbackground () const
 
Beambeam ()
 
const Beambeam () const
 
void checkNodeArgs () const
 Raises exception if a parameter value is invalid. More...
 
std::string className () const final
 Returns the class name, to be hard-coded in each leaf class that inherits from INode. More...
 
ICoordSystemcreateCoordSystem () const override
 
IDetectordetector ()
 
const IDetectordetector () const
 
bool force_polarized () const override
 Force polarized computation even in absence of sample magnetization or external fields. More...
 
const IDetectorgetDetector () const
 
const std::vector< ParameterDistribution > & getDistributions () const
 
size_t intensityMapSize () const override
 Returns the total number of the intensity values in the simulation result. More...
 
void maskAll ()
 Put the mask for all detector channels (i.e. exclude whole detector from the analysis) More...
 
std::vector< const INode * > nodeChildren () const override
 Returns all children. More...
 
std::vector< const INode * > nodeOffspring () const
 Returns all descendants. More...
 
SimulationOptionsoptions ()
 
const SimulationOptionsoptions () const
 
virtual std::vector< ParaMetaparDefs () const
 Returns the parameter definitions, to be hard-coded in each leaf class. More...
 
const MultiLayersample () const
 
void setBackground (const IBackground &bg)
 
void setBeamParameters (double lambda, int nbins, double alpha_i_min, double alpha_i_max, const IFootprintFactor *beam_shape=nullptr)
 Sets beam parameters with alpha_i of the beam defined in the range. More...
 
void setRegionOfInterest (double xlow, double ylow, double xup, double yup)
 Sets rectangular region of interest with lower left and upper right corners defined. More...
 
void setTerminalProgressMonitor ()
 Initializes a progress monitor that prints to stdout. More...
 
void setZSpan (size_t n_bins, double z_min, double z_max)
 Set z positions for intensity calculations. Negative z's correspond to the area under sample surface. The more negative z is, the deeper layer corresponds to it. More...
 
SimulationResult simulate ()
 Run a simulation, and return the result. More...
 
void subscribe (const std::function< bool(size_t)> &inform)
 
std::string unitOfParameter (ParameterDistribution::WhichParameter which) const
 
const IAxiszAxis () const
 Returns a pointer to z-position axis. More...
 

Protected Member Functions

virtual void initCoordSystem ()
 
void initDistributionHandler () override
 init callbacks for setting the parameter values More...
 
ProgressHandlerprogress ()
 
virtual void updateIntensityMap ()
 

Protected Attributes

std::unique_ptr< Beamm_beam
 
std::unique_ptr< IDetectorm_detector
 
DistributionHandler m_distribution_handler
 
std::vector< double > m_P
 

Private Member Functions

void addBackgroundIntensity (size_t start_ind, size_t n_elements) override
 
void addDataToCache (double weight) override
 
void checkCache () const
 
std::unique_ptr< IComputationcreateComputation (const reSample &re_sample, size_t start, size_t n_elements) override
 Generate a single threaded computation for a given range of simulation elements. More...
 
std::unique_ptr< DatafieldcreateIntensityData () const
 Creates intensity data from simulation elements. More...
 
std::vector< DepthProbeElementgenerateElements (const Beam &beam)
 Generate simulation elements for given beam. More...
 
double incidentAngle (size_t index) const
 
void initElementVector () override
 Initializes the vector of ISimulation elements. More...
 
void moveDataFromCache () override
 
void normalize (size_t start_ind, size_t n_elements) override
 Normalize the detector counts to beam intensity, to solid angle, and to exposure angle. More...
 
size_t numberOfElements () const override
 Gets the number of elements this simulation needs to calculate. More...
 
SimulationResult pack_result () override
 Sets m_result. More...
 
void prepareSimulation () override
 Put into a clean state for running a simulation. More...
 
void runSingleSimulation (const reSample &re_sample, size_t batch_start, size_t batch_size, double weight=1.0)
 Runs a single simulation with fixed parameter values. If desired, the simulation is run in several threads. More...
 
void setBeamParameters (double lambda, const IAxis &alpha_axis, const IFootprintFactor *beam_shape)
 Sets beam parameters with alpha_i of the beam defined in the range. More...
 
void validateParametrization (const ParameterDistribution &par_distr) const override
 Checks the distribution validity for simulation. More...
 
void validityCheck () const
 Checks if simulation data is ready for retrieval. More...
 

Private Attributes

std::unique_ptr< IAxism_alpha_axis
 
std::shared_ptr< IBackgroundm_background
 
std::vector< std::valarray< double > > m_cache
 
std::unique_ptr< DetectorContextm_detector_context
 
std::vector< DepthProbeElementm_eles
 
std::shared_ptr< SimulationOptionsm_options
 
std::shared_ptr< ProgressHandlerm_progress
 
std::shared_ptr< MultiLayerm_sample
 
std::unique_ptr< IAxism_z_axis
 

Constructor & Destructor Documentation

◆ DepthProbeSimulation()

DepthProbeSimulation::DepthProbeSimulation ( const MultiLayer sample)

Definition at line 38 of file DepthProbeSimulation.cpp.

40 {
41  // allow for negative inclinations in the beam of specular simulation
42  // it is required for proper averaging in the case of divergent beam
44 }
#define M_PI_2
Definition: Constants.h:45
void setInclinationLimits(const RealLimits &limits)
Definition: Beam.cpp:107
Beam & beam()
Definition: ISimulation2D.h:57
ISimulation2D(const Beam &beam, const MultiLayer &sample, const IDetector &detector)
const MultiLayer * sample() const
static RealLimits limited(double left_bound_value, double right_bound_value)
Creates an object bounded from the left and right.
Definition: RealLimits.cpp:134

References ISimulation2D::beam(), RealLimits::limited(), M_PI_2, and Beam::setInclinationLimits().

Here is the call graph for this function:

◆ ~DepthProbeSimulation()

DepthProbeSimulation::~DepthProbeSimulation ( )
overridedefault

Member Function Documentation

◆ addBackgroundIntensity()

void DepthProbeSimulation::addBackgroundIntensity ( size_t  start_ind,
size_t  n_elements 
)
overrideprivatevirtual

Implements ISimulation.

Definition at line 219 of file DepthProbeSimulation.cpp.

220 {
221  if (background())
222  throw std::runtime_error(
223  "Error: nonzero background is not supported by DepthProbeSimulation");
224 }
const IBackground * background() const
Definition: ISimulation.h:76

References ISimulation::background().

Here is the call graph for this function:

◆ addDataToCache()

void DepthProbeSimulation::addDataToCache ( double  weight)
overrideprivatevirtual

Implements ISimulation.

Definition at line 226 of file DepthProbeSimulation.cpp.

227 {
228  checkCache();
229  for (size_t i = 0, size = m_eles.size(); i < size; ++i)
230  m_cache[i] += m_eles[i].getIntensities() * weight;
231 }
std::vector< DepthProbeElement > m_eles
std::vector< std::valarray< double > > m_cache

References checkCache(), m_cache, and m_eles.

Here is the call graph for this function:

◆ addMask()

void ISimulation2D::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 56 of file ISimulation2D.cpp.

57 {
58  detector().addMask(shape, mask_value);
59 }
void addMask(const IShape2D &shape, bool mask_value=true)
Adds mask of given shape to the stack of detector masks. The mask value 'true' means that the channel...
Definition: IDetector.cpp:360
IDetector & detector()
Definition: ISimulation2D.h:58

References IDetector::addMask(), and ISimulation2D::detector().

Here is the call graph for this function:

◆ addParameterDistribution() [1/2]

void ISimulation::addParameterDistribution ( const ParameterDistribution par_distr)
inherited

Definition at line 231 of file ISimulation.cpp.

232 {
233  validateParametrization(par_distr);
235 }
void addParameterDistribution(const std::string &param_name, const IDistribution1D &distribution, size_t nbr_samples, double sigma_factor=0.0, const RealLimits &limits=RealLimits())
add a sampled parameter distribution
virtual void validateParametrization(const ParameterDistribution &) const
Checks the distribution validity for simulation.
Definition: ISimulation.h:120
DistributionHandler m_distribution_handler
Definition: ISimulation.h:104

References DistributionHandler::addParameterDistribution(), ISimulation::m_distribution_handler, and ISimulation::validateParametrization().

Here is the call graph for this function:

◆ addParameterDistribution() [2/2]

void ISimulation::addParameterDistribution ( ParameterDistribution::WhichParameter  whichParameter,
const IDistribution1D distribution,
size_t  nbr_samples,
double  sigma_factor = 0.0,
const RealLimits limits = RealLimits() 
)
inherited

Definition at line 222 of file ISimulation.cpp.

225 {
226  ParameterDistribution par_distr(whichParameter, distribution, nbr_samples, sigma_factor,
227  limits);
228  addParameterDistribution(par_distr);
229 }
void addParameterDistribution(ParameterDistribution::WhichParameter whichParameter, const IDistribution1D &distribution, size_t nbr_samples, double sigma_factor=0.0, const RealLimits &limits=RealLimits())
A parametric distribution function, for use with any model parameter.

◆ alphaAxis()

const IAxis * DepthProbeSimulation::alphaAxis ( ) const

Returns a pointer to incident angle axis.

Definition at line 75 of file DepthProbeSimulation.cpp.

76 {
77  if (!m_alpha_axis)
78  throw std::runtime_error("Error in DepthProbeSimulation::alphaAxis: incident angle axis "
79  "was not initialized.");
80  return m_alpha_axis.get();
81 }
std::unique_ptr< IAxis > m_alpha_axis

References m_alpha_axis.

Referenced by createIntensityData(), generateElements(), numberOfElements(), and validityCheck().

◆ background()

const IBackground* ISimulation::background ( ) const
inlineinherited

Definition at line 76 of file ISimulation.h.

76 { return m_background.get(); }
std::shared_ptr< IBackground > m_background
Definition: ISimulation.h:137

References ISimulation::m_background.

Referenced by addBackgroundIntensity(), ISimulation2D::addBackgroundIntensity(), and SpecularSimulation::addBackgroundIntensity().

◆ beam() [1/2]

◆ beam() [2/2]

const Beam& ISimulation2D::beam ( ) const
inlineinherited

Definition at line 62 of file ISimulation2D.h.

62 { return *m_beam; }

References ISimulation2D::m_beam.

◆ checkCache()

void DepthProbeSimulation::checkCache ( ) const
private

Definition at line 186 of file DepthProbeSimulation.cpp.

187 {
188  if (m_eles.size() != m_cache.size())
189  throw std::runtime_error("Error in DepthProbeSimulation: the sizes of simulation element "
190  "vector and of its cache are different");
191 }

References m_cache, and m_eles.

Referenced by addDataToCache(), and moveDataFromCache().

◆ checkNodeArgs()

void INode::checkNodeArgs ( ) const
inherited

Raises exception if a parameter value is invalid.

Definition at line 27 of file INode.cpp.

28 {
29  size_t nP = m_P.size();
30  if (parDefs().size() != nP) {
31  std::cerr << "BUG in class " << className() << std::endl;
32  std::cerr << "#m_P = " << nP << std::endl;
33  std::cerr << "#PDf = " << parDefs().size() << std::endl;
34  for (const ParaMeta& pm : parDefs())
35  std::cerr << " PDf: " << pm.name << std::endl;
36  ASSERT(0);
37  }
38  ASSERT(parDefs().size() == nP);
39  for (size_t i = 0; i < nP; ++i) {
40  const ParaMeta pm = parDefs()[i];
41 
43  if (pm.vMin == -INF) {
44  ASSERT(pm.vMax == +INF);
45  // nothing to do
46  } else if (pm.vMax == +INF) {
47  ASSERT(pm.vMin == 0);
48  limits = RealLimits::nonnegative();
49  } else {
50  limits = RealLimits::limited(pm.vMin, pm.vMax);
51  }
52  limits.check(pm.name, m_P[i]);
53  }
54 }
#define ASSERT(condition)
Definition: Assert.h:45
const double INF
Definition: INode.h:26
virtual std::vector< ParaMeta > parDefs() const
Returns the parameter definitions, to be hard-coded in each leaf class.
Definition: INode.h:51
std::vector< double > m_P
Definition: INode.h:63
virtual std::string className() const =0
Returns the class name, to be hard-coded in each leaf class that inherits from INode.
Limits for a real fit parameter.
Definition: RealLimits.h:24
static RealLimits limitless()
Creates an object without bounds (default)
Definition: RealLimits.cpp:139
void check(const std::string &name, double value) const
Throws if value is outside limits. Parameter 'name' is for exception message.
Definition: RealLimits.cpp:170
static RealLimits nonnegative()
Creates an object which can have only positive values with 0. included.
Definition: RealLimits.cpp:124
Metadata of one model parameter.
Definition: INode.h:29
double vMin
Definition: INode.h:33
double vMax
Definition: INode.h:34
std::string name
Definition: INode.h:30

References ASSERT, RealLimits::check(), INode::className(), INF, RealLimits::limited(), RealLimits::limitless(), INode::m_P, ParaMeta::name, RealLimits::nonnegative(), INode::parDefs(), ParaMeta::vMax, and ParaMeta::vMin.

Referenced by BarGauss::BarGauss(), BarLorentz::BarLorentz(), Bipyramid4::Bipyramid4(), Box::Box(), CantellatedCube::CantellatedCube(), Cone::Cone(), ConstantBackground::ConstantBackground(), CosineRippleBox::CosineRippleBox(), CosineRippleGauss::CosineRippleGauss(), CosineRippleLorentz::CosineRippleLorentz(), Cylinder::Cylinder(), DistributionCosine::DistributionCosine(), DistributionGate::DistributionGate(), DistributionGaussian::DistributionGaussian(), DistributionLogNormal::DistributionLogNormal(), DistributionLorentz::DistributionLorentz(), DistributionTrapezoid::DistributionTrapezoid(), Dodecahedron::Dodecahedron(), EllipsoidalCylinder::EllipsoidalCylinder(), FootprintGauss::FootprintGauss(), FootprintSquare::FootprintSquare(), FuzzySphere::FuzzySphere(), GaussSphere::GaussSphere(), HemiEllipsoid::HemiEllipsoid(), HollowSphere::HollowSphere(), HorizontalCylinder::HorizontalCylinder(), Icosahedron::Icosahedron(), LongBoxGauss::LongBoxGauss(), LongBoxLorentz::LongBoxLorentz(), PlatonicOctahedron::PlatonicOctahedron(), PlatonicTetrahedron::PlatonicTetrahedron(), Prism3::Prism3(), Prism6::Prism6(), Profile1DCauchy::Profile1DCauchy(), Profile1DCosine::Profile1DCosine(), Profile1DGate::Profile1DGate(), Profile1DGauss::Profile1DGauss(), Profile1DTriangle::Profile1DTriangle(), Profile1DVoigt::Profile1DVoigt(), Profile2DCauchy::Profile2DCauchy(), Profile2DCone::Profile2DCone(), Profile2DGate::Profile2DGate(), Profile2DGauss::Profile2DGauss(), Profile2DVoigt::Profile2DVoigt(), Pyramid2::Pyramid2(), Pyramid3::Pyramid3(), Pyramid4::Pyramid4(), Pyramid6::Pyramid6(), RotationEuler::RotationEuler(), RotationX::RotationX(), RotationY::RotationY(), RotationZ::RotationZ(), SawtoothRippleBox::SawtoothRippleBox(), SawtoothRippleGauss::SawtoothRippleGauss(), SawtoothRippleLorentz::SawtoothRippleLorentz(), Sphere::Sphere(), Spheroid::Spheroid(), TruncatedCube::TruncatedCube(), TruncatedSphere::TruncatedSphere(), and TruncatedSpheroid::TruncatedSpheroid().

Here is the call graph for this function:

◆ className()

std::string DepthProbeSimulation::className ( ) const
inlinefinalvirtual

Returns the class name, to be hard-coded in each leaf class that inherits from INode.

Implements INode.

Definition at line 42 of file DepthProbeSimulation.h.

42 { return "DepthProbeSimulation"; }

◆ createComputation()

std::unique_ptr< IComputation > DepthProbeSimulation::createComputation ( const reSample re_sample,
size_t  start,
size_t  n_elements 
)
overrideprivatevirtual

Generate a single threaded computation for a given range of simulation elements.

Parameters
re_samplePreprocessed version of our sample
startIndex of the first element to include into computation
n_elementsNumber of elements to process

Implements ISimulation.

Definition at line 164 of file DepthProbeSimulation.cpp.

165 {
166  ASSERT(start < m_eles.size() && start + n_elements <= m_eles.size());
167  const auto& begin = m_eles.begin() + static_cast<long>(start);
168  return std::make_unique<DepthProbeComputation>(re_sample, options(), progress(), begin,
169  begin + static_cast<long>(n_elements));
170 }
ProgressHandler & progress()
const SimulationOptions & options() const

References ASSERT, m_eles, ISimulation::options(), and ISimulation::progress().

Here is the call graph for this function:

◆ createCoordSystem()

ICoordSystem * DepthProbeSimulation::createCoordSystem ( ) const
overridevirtual

Implements ISimulation.

Definition at line 100 of file DepthProbeSimulation.cpp.

101 {
102  OwningVector<IAxis> axes({m_alpha_axis->clone(), m_z_axis->clone()});
103  return new DepthProbeCoordinates(axes, beam().direction(), beam().wavelength());
104 }
DepthProbeCoordinates class handles the unit translations for depth probe simulations Its default uni...
std::unique_ptr< IAxis > m_z_axis

References ISimulation2D::beam(), m_alpha_axis, and m_z_axis.

Referenced by pack_result().

Here is the call graph for this function:

◆ createIntensityData()

std::unique_ptr< Datafield > DepthProbeSimulation::createIntensityData ( ) const
private

Creates intensity data from simulation elements.

Definition at line 247 of file DepthProbeSimulation.cpp.

248 {
249  auto frame = new Frame({alphaAxis()->clone(), zAxis()->clone()});
250 
251  std::vector<double> out(frame->size());
252  size_t iout = 0;
253  for (size_t i = 0, size = m_eles.size(); i < size; ++i)
254  for (const double v: m_eles[i].getIntensities())
255  out[iout++] = v;
256 
257  return std::unique_ptr<Datafield>(new Datafield(frame, out));
258 }
Stores radiation power per bin.
Definition: Datafield.h:30
const IAxis * alphaAxis() const
Returns a pointer to incident angle axis.
const IAxis * zAxis() const
Returns a pointer to z-position axis.
Holds one or two axes.
Definition: Frame.h:27
virtual IAxis * clone() const =0

References alphaAxis(), IAxis::clone(), m_eles, and zAxis().

Referenced by pack_result().

Here is the call graph for this function:

◆ detector() [1/2]

◆ detector() [2/2]

const IDetector& ISimulation2D::detector ( ) const
inlineinherited

Definition at line 63 of file ISimulation2D.h.

63 { return *m_detector; }

References ISimulation2D::m_detector.

◆ force_polarized()

bool DepthProbeSimulation::force_polarized ( ) const
inlineoverridevirtual

Force polarized computation even in absence of sample magnetization or external fields.

Implements ISimulation.

Definition at line 61 of file DepthProbeSimulation.h.

61 { return false; }

◆ generateElements()

std::vector< DepthProbeElement > DepthProbeSimulation::generateElements ( const Beam beam)
private

Generate simulation elements for given beam.

Definition at line 145 of file DepthProbeSimulation.cpp.

146 {
147  std::vector<DepthProbeElement> result;
148 
149  const double wavelength = beam.wavelength();
150  const double angle_shift = beam.direction().alpha();
151 
152  const size_t axis_size = alphaAxis()->size();
153  result.reserve(axis_size);
154  for (size_t i = 0; i < axis_size; ++i) {
155  double result_angle = incidentAngle(i) + angle_shift;
156  result.emplace_back(wavelength, -result_angle, zAxis());
157  if (!alpha_limits.isInRange(result_angle))
158  result.back().setCalculationFlag(false); // false = exclude from calculations
159  }
160  return result;
161 }
Direction direction() const
Definition: Beam.h:46
double wavelength() const
Definition: Beam.h:44
double incidentAngle(size_t index) const
double alpha() const
Definition: Direction.h:36
virtual size_t size() const =0
Returns the number of bins.

References Direction::alpha(), alphaAxis(), ISimulation2D::beam(), Beam::direction(), incidentAngle(), IAxis::size(), Beam::wavelength(), and zAxis().

Referenced by initElementVector().

Here is the call graph for this function:

◆ getDetector()

const IDetector* ISimulation2D::getDetector ( ) const
inlineinherited

Definition at line 64 of file ISimulation2D.h.

64 { return m_detector.get(); }

References ISimulation2D::m_detector.

Referenced by ScatteringSimulation::createCoordSystem().

◆ getDistributions()

const std::vector< ParameterDistribution > & ISimulation::getDistributions ( ) const
inherited

Definition at line 237 of file ISimulation.cpp.

238 {
240 }
const std::vector< ParameterDistribution > & getDistributions() const

References DistributionHandler::getDistributions(), and ISimulation::m_distribution_handler.

Here is the call graph for this function:

◆ incidentAngle()

double DepthProbeSimulation::incidentAngle ( size_t  index) const
private

Definition at line 242 of file DepthProbeSimulation.cpp.

243 {
244  return m_alpha_axis->bin(index).center();
245 }

References m_alpha_axis.

Referenced by generateElements().

◆ initCoordSystem()

virtual void ISimulation2D::initCoordSystem ( )
inlineprotectedvirtualinherited

Definition at line 70 of file ISimulation2D.h.

70 {}

◆ initDistributionHandler()

void ISimulation2D::initDistributionHandler ( )
overrideprotectedvirtualinherited

init callbacks for setting the parameter values

Reimplemented from ISimulation.

Definition at line 159 of file ISimulation2D.cpp.

160 {
161  for (const auto& distribution : m_distribution_handler.getDistributions()) {
162  const auto* dc = const_cast<ParameterDistribution*>(&distribution);
163 
164  switch (dc->whichParameter()) {
167  dc, [&](double d) { beam().setAzimuthalAngleGuarded(d); });
168  break;
171  dc, [&](double d) { beam().setInclinationAngleGuarded(d); });
172  break;
175  dc, [&](double d) { beam().setWavelengthGuarded(d); });
176  break;
177  default:
178  ASSERT(0);
179  }
180  }
181 }
void setAzimuthalAngleGuarded(double value)
Check for limits, set the value if within limits. Throws if limits are violated.
Definition: Beam.cpp:123
void setInclinationAngleGuarded(double value)
Check for limits, set the value if within limits. Throws if limits are violated.
Definition: Beam.cpp:117
void setWavelengthGuarded(double value)
Check for limits, set the value if within limits. Throws if limits are violated.
Definition: Beam.cpp:129
void defineCallbackForDistribution(const ParameterDistribution *distribution, std::function< void(double)> fn)

References ASSERT, ISimulation2D::beam(), ParameterDistribution::BeamAzimuthalAngle, ParameterDistribution::BeamInclinationAngle, ParameterDistribution::BeamWavelength, DistributionHandler::defineCallbackForDistribution(), DistributionHandler::getDistributions(), ISimulation::m_distribution_handler, Beam::setAzimuthalAngleGuarded(), Beam::setInclinationAngleGuarded(), and Beam::setWavelengthGuarded().

Here is the call graph for this function:

◆ initElementVector()

void DepthProbeSimulation::initElementVector ( )
overrideprivatevirtual

Initializes the vector of ISimulation elements.

Implements ISimulation.

Definition at line 136 of file DepthProbeSimulation.cpp.

137 {
139 
140  if (!m_cache.empty())
141  return;
142  m_cache.resize(m_eles.size(), std::valarray<double>(0.0, zAxis()->size()));
143 }
std::vector< DepthProbeElement > generateElements(const Beam &beam)
Generate simulation elements for given beam.

References ISimulation2D::beam(), generateElements(), m_cache, m_eles, and zAxis().

Here is the call graph for this function:

◆ intensityMapSize()

size_t DepthProbeSimulation::intensityMapSize ( ) const
overridevirtual

Returns the total number of the intensity values in the simulation result.

Implements ISimulation.

Definition at line 91 of file DepthProbeSimulation.cpp.

92 {
93  if (!m_z_axis || !m_alpha_axis)
94  throw std::runtime_error("Error in DepthProbeSimulation::intensityMapSize: attempt to "
95  "access non-initialized data.");
96  return m_z_axis->size() * m_alpha_axis->size();
97 }

References m_alpha_axis, and m_z_axis.

◆ maskAll()

void ISimulation2D::maskAll ( )
inherited

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

Definition at line 61 of file ISimulation2D.cpp.

62 {
63  detector().maskAll();
64 }
void maskAll()
Put the mask for all detector channels (i.e. exclude whole detector from the analysis)
Definition: IDetector.cpp:365

References ISimulation2D::detector(), and IDetector::maskAll().

Here is the call graph for this function:

◆ moveDataFromCache()

void DepthProbeSimulation::moveDataFromCache ( )
overrideprivatevirtual

Implements ISimulation.

Definition at line 233 of file DepthProbeSimulation.cpp.

234 {
235  checkCache();
236  for (size_t i = 0, size = m_eles.size(); i < size; ++i)
237  m_eles[i].setIntensities(std::move(m_cache[i]));
238  m_cache.clear();
239  m_cache.shrink_to_fit();
240 }

References checkCache(), m_cache, and m_eles.

Here is the call graph for this function:

◆ nodeChildren()

std::vector< const INode * > ISimulation2D::nodeChildren ( ) const
overridevirtualinherited

Returns all children.

Reimplemented from ISimulation.

Definition at line 42 of file ISimulation2D.cpp.

43 {
44  std::vector<const INode*> result = ISimulation::nodeChildren();
45  result.push_back(m_beam.get());
46  if (m_detector)
47  result.push_back(m_detector.get());
48  return result;
49 }
std::vector< const INode * > nodeChildren() const override
Returns all children.

References ISimulation2D::m_beam, ISimulation2D::m_detector, and ISimulation::nodeChildren().

Here is the call graph for this function:

◆ nodeOffspring()

std::vector< const INode * > INode::nodeOffspring ( ) const
inherited

Returns all descendants.

Definition at line 61 of file INode.cpp.

62 {
63  std::vector<const INode*> result;
64  result.push_back(this);
65  for (const auto* child : nodeChildren()) {
66  for (const auto* p : child->nodeOffspring())
67  result.push_back(p);
68  }
69  return result;
70 }
virtual std::vector< const INode * > nodeChildren() const
Returns all children.
Definition: INode.cpp:56

References INode::nodeChildren().

Here is the call graph for this function:

◆ normalize()

void DepthProbeSimulation::normalize ( size_t  start_ind,
size_t  n_elements 
)
overrideprivatevirtual

Normalize the detector counts to beam intensity, to solid angle, and to exposure angle.

Parameters
start_indIndex of the first element to operate on
n_elementsNumber of elements to process

Implements ISimulation.

Definition at line 205 of file DepthProbeSimulation.cpp.

206 {
207  const double beam_intensity = beam().intensity();
208  for (size_t i = start_ind, stop_point = start_ind + n_elements; i < stop_point; ++i) {
209  auto& element = m_eles[i];
210  const double alpha_i = -element.alphaI();
211  const auto* footprint = beam().footprintFactor();
212  double intensity_factor = beam_intensity;
213  if (footprint != nullptr)
214  intensity_factor = intensity_factor * footprint->calculate(alpha_i);
215  element.setIntensities(element.getIntensities() * intensity_factor);
216  }
217 }
double intensity() const
Returns the beam intensity in neutrons/sec.
Definition: Beam.h:43
const IFootprintFactor * footprintFactor() const
Returns footprint factor.
Definition: Beam.cpp:112
virtual double calculate(double alpha) const =0
Calculate footprint correction coefficient from the beam incident angle alpha.

References ISimulation2D::beam(), IFootprintFactor::calculate(), Beam::footprintFactor(), Beam::intensity(), and m_eles.

Here is the call graph for this function:

◆ numberOfElements()

size_t DepthProbeSimulation::numberOfElements ( ) const
overrideprivatevirtual

Gets the number of elements this simulation needs to calculate.

Implements ISimulation.

Definition at line 48 of file DepthProbeSimulation.cpp.

49 {
50  return alphaAxis()->size();
51 }

References alphaAxis(), and IAxis::size().

Here is the call graph for this function:

◆ options() [1/2]

SimulationOptions & ISimulation::options ( )
inherited

Definition at line 123 of file ISimulation.cpp.

124 {
125  ASSERT(m_options);
126  return *m_options;
127 }
std::shared_ptr< SimulationOptions > m_options
Definition: ISimulation.h:135

References ASSERT, and ISimulation::m_options.

◆ options() [2/2]

const SimulationOptions & ISimulation::options ( ) const
inherited

◆ pack_result()

SimulationResult DepthProbeSimulation::pack_result ( )
overrideprivatevirtual

Sets m_result.

Implements ISimulation.

Definition at line 53 of file DepthProbeSimulation.cpp.

54 {
55  validityCheck();
56  auto data = createIntensityData();
57  return {*data, createCoordSystem()};
58 }
std::unique_ptr< Datafield > createIntensityData() const
Creates intensity data from simulation elements.
ICoordSystem * createCoordSystem() const override
void validityCheck() const
Checks if simulation data is ready for retrieval.

References createCoordSystem(), createIntensityData(), and validityCheck().

Here is the call graph for this function:

◆ parDefs()

virtual std::vector<ParaMeta> INode::parDefs ( ) const
inlinevirtualinherited

Returns the parameter definitions, to be hard-coded in each leaf class.

Reimplemented in ConstantBackground, GaussSphere, FuzzySphere, RotationEuler, RotationZ, RotationY, RotationX, Crystal, Layer, HexagonalLattice2D, SquareLattice2D, BasicLattice2D, LayerRoughness, TruncatedSpheroid, TruncatedSphere, TruncatedCube, Spheroid, Sphere, SawtoothRippleLorentz, SawtoothRippleGauss, SawtoothRippleBox, Pyramid6, Pyramid4, Pyramid3, Pyramid2, Prism6, Prism3, PlatonicTetrahedron, PlatonicOctahedron, LongBoxLorentz, LongBoxGauss, Icosahedron, HorizontalCylinder, HollowSphere, HemiEllipsoid, EllipsoidalCylinder, Dodecahedron, Cylinder, CosineRippleLorentz, CosineRippleGauss, CosineRippleBox, Cone, CantellatedCube, Box, Bipyramid4, BarLorentz, BarGauss, Profile2DVoigt, Profile2DCone, Profile2DGate, Profile2DGauss, Profile2DCauchy, Profile1DVoigt, Profile1DCosine, Profile1DTriangle, Profile1DGate, Profile1DGauss, Profile1DCauchy, MisesGaussPeakShape, MisesFisherGaussPeakShape, LorentzFisherPeakShape, GaussFisherPeakShape, IsotropicLorentzPeakShape, IsotropicGaussPeakShape, ParticleLayout, InterferenceTwin, InterferenceRadialParaCrystal, InterferenceHardDisk, Interference2DSuperLattice, Interference2DParaCrystal, Interference1DLattice, DistributionTrapezoid, DistributionCosine, DistributionLogNormal, DistributionGaussian, DistributionLorentz, DistributionGate, ResolutionFunction2DGaussian, PolFilter, FootprintSquare, and FootprintGauss.

Definition at line 51 of file INode.h.

51 { return {}; }

Referenced by INode::checkNodeArgs(), and IFormFactor::pythonConstructor().

◆ prepareSimulation()

void DepthProbeSimulation::prepareSimulation ( )
inlineoverrideprivatevirtual

Put into a clean state for running a simulation.

Implements ISimulation.

Definition at line 67 of file DepthProbeSimulation.h.

67 {}

◆ progress()

ProgressHandler & ISimulation::progress ( )
protectedinherited

Definition at line 129 of file ISimulation.cpp.

130 {
132  return *m_progress;
133 }
std::shared_ptr< ProgressHandler > m_progress
Definition: ISimulation.h:136

References ASSERT, and ISimulation::m_progress.

Referenced by createComputation(), ISimulation2D::createComputation(), and SpecularSimulation::createComputation().

◆ runSingleSimulation()

void ISimulation::runSingleSimulation ( const reSample re_sample,
size_t  batch_start,
size_t  batch_size,
double  weight = 1.0 
)
privateinherited

Runs a single simulation with fixed parameter values. If desired, the simulation is run in several threads.

Definition at line 245 of file ISimulation.cpp.

247 {
249 
250  const size_t n_threads = m_options->getNumberOfThreads();
251  ASSERT(n_threads > 0);
252 
253  std::vector<std::unique_ptr<IComputation>> computations;
254 
255  for (size_t i_thread = 0; i_thread < n_threads; ++i_thread) {
256  const size_t thread_start = batch_start + startIndex(n_threads, i_thread, batch_size);
257  const size_t thread_size = batchSize(n_threads, i_thread, batch_size);
258  if (thread_size == 0)
259  break;
260  computations.emplace_back(createComputation(re_sample, thread_start, thread_size));
261  }
262  runComputations(computations);
263 
264  normalize(batch_start, batch_size);
265  addBackgroundIntensity(batch_start, batch_size);
266  addDataToCache(weight);
267 }
virtual void initElementVector()=0
Initializes the vector of ISimulation elements.
virtual void addDataToCache(double weight)=0
virtual void addBackgroundIntensity(size_t start_ind, size_t n_elements)=0
virtual std::unique_ptr< IComputation > createComputation(const reSample &re_sample, size_t start, size_t n_elements)=0
Generate a single threaded computation for a given range of simulation elements.
virtual void normalize(size_t start_ind, size_t n_elements)=0
Normalize the detector counts to beam intensity, to solid angle, and to exposure angle.

References ISimulation::addBackgroundIntensity(), ISimulation::addDataToCache(), ASSERT, ISimulation::createComputation(), ISimulation::initElementVector(), ISimulation::m_options, and ISimulation::normalize().

Referenced by ISimulation::simulate().

Here is the call graph for this function:

◆ sample()

const MultiLayer * ISimulation::sample ( ) const
inherited

Definition at line 191 of file ISimulation.cpp.

192 {
193  return m_sample.get();
194 }
std::shared_ptr< MultiLayer > m_sample
Definition: ISimulation.h:133

References ISimulation::m_sample.

Referenced by ISimulation::simulate(), and validityCheck().

◆ setBackground()

void ISimulation::setBackground ( const IBackground bg)
inherited

Definition at line 196 of file ISimulation.cpp.

197 {
198  m_background.reset(bg.clone());
199 }
IBackground * clone() const override=0

References IBackground::clone(), and ISimulation::m_background.

Here is the call graph for this function:

◆ setBeamParameters() [1/2]

void DepthProbeSimulation::setBeamParameters ( double  lambda,
const IAxis alpha_axis,
const IFootprintFactor beam_shape 
)
private

Sets beam parameters with alpha_i of the beam defined in the range.

Definition at line 107 of file DepthProbeSimulation.cpp.

109 {
110  if (lambda <= 0.0)
111  throw std::runtime_error(
112  "Error in DepthProbeSimulation::setBeamParameters: wavelength must be positive.");
113  if (alpha_axis.min() < 0.0)
114  throw std::runtime_error(
115  "Error in DepthProbeSimulation::setBeamParameters: minimum value on "
116  "angle axis is negative.");
117  if (alpha_axis.min() >= alpha_axis.max())
118  throw std::runtime_error(
119  "Error in DepthProbeSimulation::setBeamParameters: alpha bounds are not sorted.");
120  if (alpha_axis.size() == 0)
121  throw std::runtime_error(
122  "Error in DepthProbeSimulation::setBeamParameters: angle axis is empty");
123 
124  m_alpha_axis.reset(alpha_axis.clone());
125 
126  // beam is initialized with zero-valued angles
127  // Zero-valued incident alpha is required for proper
128  // taking into account beam resolution effects
129  beam().setWavelength(lambda);
130  beam().setDirection({zero_alpha_i, zero_phi_i});
131 
132  if (beam_shape)
133  beam().setFootprintFactor(*beam_shape);
134 }
void setDirection(const Direction &direction)
Definition: Beam.cpp:92
void setFootprintFactor(const IFootprintFactor &shape_factor)
Sets footprint factor to the beam.
Definition: Beam.cpp:135
void setWavelength(double wavelength)
Definition: Beam.cpp:85
virtual double max() const =0
Returns value of last point of axis.
virtual double min() const =0
Returns value of first point of axis.

References ISimulation2D::beam(), IAxis::clone(), m_alpha_axis, IAxis::max(), IAxis::min(), Beam::setDirection(), Beam::setFootprintFactor(), Beam::setWavelength(), and IAxis::size().

Here is the call graph for this function:

◆ setBeamParameters() [2/2]

void DepthProbeSimulation::setBeamParameters ( double  lambda,
int  nbins,
double  alpha_i_min,
double  alpha_i_max,
const IFootprintFactor beam_shape = nullptr 
)

Sets beam parameters with alpha_i of the beam defined in the range.

Definition at line 60 of file DepthProbeSimulation.cpp.

62 {
63  FixedBinAxis axis("alpha_i", static_cast<size_t>(nbins), alpha_i_min, alpha_i_max);
64  setBeamParameters(lambda, axis, beam_shape);
65 }
void setBeamParameters(double lambda, int nbins, double alpha_i_min, double alpha_i_max, const IFootprintFactor *beam_shape=nullptr)
Sets beam parameters with alpha_i of the beam defined in the range.
Axis with fixed bin size.
Definition: FixedBinAxis.h:23

◆ setRegionOfInterest()

void ISimulation2D::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 66 of file ISimulation2D.cpp.

67 {
68  detector().setRegionOfInterest(xlow, ylow, xup, yup);
69 }
void setRegionOfInterest(double xlow, double ylow, double xup, double yup)
Sets rectangular region of interest with lower left and upper right corners defined.
Definition: IDetector.cpp:337

References ISimulation2D::detector(), and IDetector::setRegionOfInterest().

Here is the call graph for this function:

◆ setTerminalProgressMonitor()

void ISimulation::setTerminalProgressMonitor ( )
inherited

Initializes a progress monitor that prints to stdout.

Definition at line 142 of file ISimulation.cpp.

143 {
144  m_progress->subscribe([](size_t percentage_done) -> bool {
145  if (percentage_done < 100)
146  std::cout << std::setprecision(2) << "\r... " << percentage_done << "%" << std::flush;
147  else // wipe out
148  std::cout << "\r... 100%\n";
149  return true;
150  });
151 }

References ISimulation::m_progress.

◆ setZSpan()

void DepthProbeSimulation::setZSpan ( size_t  n_bins,
double  z_min,
double  z_max 
)

Set z positions for intensity calculations. Negative z's correspond to the area under sample surface. The more negative z is, the deeper layer corresponds to it.

Definition at line 67 of file DepthProbeSimulation.cpp.

68 {
69  if (z_max <= z_min)
70  throw std::runtime_error("Error in DepthProbeSimulation::setZSpan: maximum on-axis value "
71  "is less or equal to the minimum one");
72  m_z_axis = std::make_unique<FixedBinAxis>("z", n_bins, z_min, z_max);
73 }

References m_z_axis.

◆ simulate()

SimulationResult ISimulation::simulate ( )
inherited

Run a simulation, and return the result.

Runs simulation with possible averaging over parameter distributions; returns result.

Definition at line 154 of file ISimulation.cpp.

155 {
157  throw std::runtime_error("Sample contains incompatible material definitions");
158  gsl_set_error_handler_off();
159 
161 
162  const auto re_sample = reSample::make(*sample(), options(), force_polarized());
163 
164  const size_t total_size = numberOfElements();
165  size_t param_combinations = m_distribution_handler.getTotalNumberOfSamples();
166 
167  m_progress->reset();
168  m_progress->setExpectedNTicks(param_combinations * total_size);
169 
170  // restrict calculation to current batch
171  const size_t n_batches = m_options->getNumberOfBatches();
172  const size_t current_batch = m_options->getCurrentBatch();
173 
174  const size_t batch_start = startIndex(n_batches, current_batch, total_size);
175  const size_t batch_size = batchSize(n_batches, current_batch, total_size);
176  ASSERT(batch_size);
177 
178  if (param_combinations == 1)
179  runSingleSimulation(re_sample, batch_start, batch_size, 1.);
180  else {
182  for (size_t index = 0; index < param_combinations; ++index) {
183  double weight = m_distribution_handler.setParameterValues(index);
184  runSingleSimulation(re_sample, batch_start, batch_size, weight);
185  }
186  }
188  return pack_result();
189 }
double setParameterValues(size_t index)
set the parameter values of the simulation object to a specific combination of values,...
size_t getTotalNumberOfSamples() const
get the total number of parameter value combinations (product of the individual sizes of each paramet...
virtual size_t numberOfElements() const =0
Gets the number of elements this simulation needs to calculate.
virtual void prepareSimulation()=0
Put into a clean state for running a simulation.
virtual SimulationResult pack_result()=0
Sets m_result.
virtual void moveDataFromCache()=0
void runSingleSimulation(const reSample &re_sample, size_t batch_start, size_t batch_size, double weight=1.0)
Runs a single simulation with fixed parameter values. If desired, the simulation is run in several th...
virtual bool force_polarized() const =0
Force polarized computation even in absence of sample magnetization or external fields.
virtual void initDistributionHandler()
Definition: ISimulation.h:106
static reSample make(const MultiLayer &sample, const SimulationOptions &options, bool forcePolarized=false)
Factory method that wraps the private constructor.
Definition: ReSample.cpp:309
bool ContainsCompatibleMaterials(const MultiLayer &sample)
Returns true if the sample contains non-default materials of one type only.

References ASSERT, SampleUtils::Multilayer::ContainsCompatibleMaterials(), ISimulation::force_polarized(), DistributionHandler::getTotalNumberOfSamples(), ISimulation::initDistributionHandler(), ISimulation::m_distribution_handler, ISimulation::m_options, ISimulation::m_progress, ISimulation::m_sample, reSample::make(), ISimulation::moveDataFromCache(), ISimulation::numberOfElements(), ISimulation::options(), ISimulation::pack_result(), ISimulation::prepareSimulation(), ISimulation::runSingleSimulation(), ISimulation::sample(), and DistributionHandler::setParameterValues().

Here is the call graph for this function:

◆ subscribe()

void ISimulation::subscribe ( const std::function< bool(size_t)> &  inform)
inherited

Definition at line 135 of file ISimulation.cpp.

136 {
138  m_progress->subscribe(inform);
139 }

References ASSERT, and ISimulation::m_progress.

◆ unitOfParameter()

std::string ISimulation::unitOfParameter ( ParameterDistribution::WhichParameter  which) const
inherited

Definition at line 209 of file ISimulation.cpp.

210 {
211  switch (which) {
213  return "nm";
216  return "rad";
217  default:
218  throw std::runtime_error("ISimulation::unitOfParameter(): missing implementation");
219  }
220 }

References ParameterDistribution::BeamAzimuthalAngle, ParameterDistribution::BeamInclinationAngle, and ParameterDistribution::BeamWavelength.

◆ updateIntensityMap()

virtual void ISimulation::updateIntensityMap ( )
inlineprotectedvirtualinherited

Reimplemented in OffspecSimulation.

Definition at line 97 of file ISimulation.h.

97 {}

◆ validateParametrization()

void DepthProbeSimulation::validateParametrization ( const ParameterDistribution par_distr) const
overrideprivatevirtual

Checks the distribution validity for simulation.

Reimplemented from ISimulation.

Definition at line 193 of file DepthProbeSimulation.cpp.

194 {
195  const bool zero_mean = par_distr.getDistribution()->mean() == 0.0;
196  if (zero_mean)
197  return;
198  /* TODO come back to this after refactoring of distribution handling
199  if (par_distr.whichParameter() == ParameterDistribution::BeamInclinationAngle)
200  throw std::runtime_error("Error in DepthProbeSimulation: parameter distribution of "
201  "beam inclination angle should have zero mean.");
202  */
203 }
virtual double mean() const =0
Returns the distribution-specific mean.
const IDistribution1D * getDistribution() const

References ParameterDistribution::getDistribution(), and IDistribution1D::mean().

Here is the call graph for this function:

◆ validityCheck()

void DepthProbeSimulation::validityCheck ( ) const
private

Checks if simulation data is ready for retrieval.

Definition at line 172 of file DepthProbeSimulation.cpp.

173 {
174  const MultiLayer* current_sample = sample();
175  if (!current_sample)
176  throw std::runtime_error(
177  "Error in DepthProbeSimulation::validityCheck: no sample found in the simulation.");
178 
179  const size_t data_size = m_eles.size();
180  if (data_size != alphaAxis()->size())
181  throw std::runtime_error(
182  "Error in DepthProbeSimulation::validityCheck: length of simulation "
183  "element vector is not equal to the number of inclination angles");
184 }
Our sample model: a stack of layers one below the other.
Definition: MultiLayer.h:43

References alphaAxis(), m_eles, and ISimulation::sample().

Referenced by pack_result().

Here is the call graph for this function:

◆ zAxis()

const IAxis * DepthProbeSimulation::zAxis ( ) const

Returns a pointer to z-position axis.

Definition at line 83 of file DepthProbeSimulation.cpp.

84 {
85  if (!m_z_axis)
86  throw std::runtime_error("Error in DepthProbeSimulation::zAxis: position axis "
87  "was not initialized.");
88  return m_z_axis.get();
89 }

References m_z_axis.

Referenced by createIntensityData(), generateElements(), and initElementVector().

Member Data Documentation

◆ m_alpha_axis

std::unique_ptr<IAxis> DepthProbeSimulation::m_alpha_axis
private

◆ m_background

std::shared_ptr<IBackground> ISimulation::m_background
privateinherited

Definition at line 137 of file ISimulation.h.

Referenced by ISimulation::background(), and ISimulation::setBackground().

◆ m_beam

std::unique_ptr<Beam> ISimulation2D::m_beam
protectedinherited

◆ m_cache

std::vector<std::valarray<double> > DepthProbeSimulation::m_cache
private

◆ m_detector

std::unique_ptr<IDetector> ISimulation2D::m_detector
protectedinherited

◆ m_detector_context

std::unique_ptr<DetectorContext> ISimulation2D::m_detector_context
privateinherited

◆ m_distribution_handler

DistributionHandler ISimulation::m_distribution_handler
protectedinherited

◆ m_eles

std::vector<DepthProbeElement> DepthProbeSimulation::m_eles
private

◆ m_options

std::shared_ptr<SimulationOptions> ISimulation::m_options
privateinherited

◆ m_P

◆ m_progress

std::shared_ptr<ProgressHandler> ISimulation::m_progress
privateinherited

◆ m_sample

std::shared_ptr<MultiLayer> ISimulation::m_sample
privateinherited

◆ m_z_axis

std::unique_ptr<IAxis> DepthProbeSimulation::m_z_axis
private

Definition at line 116 of file DepthProbeSimulation.h.

Referenced by createCoordSystem(), intensityMapSize(), setZSpan(), and zAxis().


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