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

Public Member Functions

 Simulation ()
 
virtual ~Simulation ()
 
virtual Simulationclone () const =0
 
virtual void prepareSimulation ()
 
void runSimulation ()
 
void runMPISimulation ()
 
void setInstrument (const Instrument &instrument_)
 
const Instrumentinstrument () const
 
Instrumentinstrument ()
 
void setBeamIntensity (double intensity)
 
double getBeamIntensity () const
 
void setBeamPolarization (const kvector_t bloch_vector)
 
void setDetectorResolutionFunction (const IResolutionFunction2D &resolution_function)
 
void removeDetectorResolutionFunction ()
 
void setAnalyzerProperties (const kvector_t direction, double efficiency, double total_transmission)
 
void setSample (const MultiLayer &sample)
 
const MultiLayersample () const
 
void setSampleBuilder (const std::shared_ptr< ISampleBuilder > &sample_builder)
 
void setBackground (const IBackground &bg)
 
const IBackgroundbackground () const
 
virtual size_t intensityMapSize () const =0
 
virtual SimulationResult result () const =0
 
void addParameterDistribution (const std::string &param_name, const IDistribution1D &distribution, size_t nbr_samples, double sigma_factor=0.0, const RealLimits &limits=RealLimits())
 
void addParameterDistribution (const ParameterDistribution &par_distr)
 
const DistributionHandlergetDistributionHandler () const
 
void setOptions (const SimulationOptions &options)
 
const SimulationOptionsgetOptions () const
 
SimulationOptionsgetOptions ()
 
void subscribe (ProgressHandler::Callback_t inform)
 
void setTerminalProgressMonitor ()
 
std::vector< const INode * > getChildren () const
 
SimulationResult convertData (const OutputData< double > &data, bool put_masked_areas_to_zero=true)
 
virtual void transferToCPP ()
 
virtual void accept (INodeVisitor *visitor) const =0
 
virtual std::string treeToString () const
 
void registerChild (INode *node)
 
virtual void setParent (const INode *newParent)
 
const INodeparent () const
 
INodeparent ()
 
int copyNumber (const INode *node) const
 
std::string displayName () const
 
ParameterPoolcreateParameterTree () const
 
ParameterPoolparameterPool () const
 
std::string parametersToString () const
 
RealParameterregisterParameter (const std::string &name, double *parpointer)
 
void registerVector (const std::string &base_name, kvector_t *p_vec, const std::string &units="nm")
 
void setParameterValue (const std::string &name, double value)
 
void setVectorValue (const std::string &base_name, kvector_t value)
 
RealParameterparameter (const std::string &name) const
 
virtual void onChange ()
 
void removeParameter (const std::string &name)
 
void removeVector (const std::string &base_name)
 
void setName (const std::string &name)
 
const std::string & getName () const
 

Static Public Member Functions

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

Protected Member Functions

 Simulation (const Simulation &other)
 
virtual void transferResultsToIntensityMap ()
 
virtual void initSimulationElementVector ()=0
 
virtual void updateIntensityMap ()
 
virtual size_t numberOfSimulationElements () const =0
 
const SimulationOptionsoptions () const
 
ProgressHandlerprogress ()
 

Protected Attributes

const size_t m_NP
 
std::vector< double > m_P
 

Private Member Functions

void initialize ()
 
void runSingleSimulation (size_t batch_start, size_t batch_size, double weight=1.0)
 
virtual std::unique_ptr< IComputationgenerateSingleThreadedComputation (size_t start, size_t n_elements)=0
 
virtual void validateParametrization (const ParameterDistribution &) const
 
virtual void addBackgroundIntensity (size_t start_ind, size_t n_elements)=0
 
virtual void normalize (size_t start_ind, size_t n_elements)=0
 
virtual void addDataToCache (double weight)=0
 
virtual void moveDataFromCache ()=0
 
virtual std::vector< double > rawResults () const =0
 
virtual void setRawResults (const std::vector< double > &raw_data)=0
 

Private Attributes

SimulationOptions m_options
 
ProgressHandler m_progress
 
SampleProvider m_sample_provider
 
DistributionHandler m_distribution_handler
 
Instrument m_instrument
 
std::unique_ptr< IBackgroundm_background
 
const INodem_parent {nullptr}
 
std::string m_name
 
std::unique_ptr< ParameterPoolm_pool
 

Friends

class MPISimulation
 

Detailed Description

Pure virtual base class of OffSpecularSimulation, GISASSimulation and SpecularSimulation.

Holds the common infrastructure to run a simulation: multithreading, batch processing, weighting over parameter distributions, ...

Definition at line 37 of file Simulation.h.

Constructor & Destructor Documentation

◆ Simulation() [1/2]

Simulation::Simulation ( )

Definition at line 121 of file Simulation.cpp.

122 {
123  initialize();
124 }
void initialize()
Definition: Simulation.cpp:138

References initialize().

Here is the call graph for this function:

◆ ~Simulation()

Simulation::~Simulation ( )
virtualdefault

◆ Simulation() [2/2]

Simulation::Simulation ( const Simulation other)
protected

Definition at line 126 of file Simulation.cpp.

127  : ICloneable(), INode(), m_options(other.m_options), m_progress(other.m_progress),
130 {
131  if (other.m_background)
132  setBackground(*other.m_background);
133  initialize();
134 }
ICloneable()=default
INode()
Definition: INode.h:51
SimulationOptions m_options
Definition: Simulation.h:152
Instrument m_instrument
Definition: Simulation.h:156
void setBackground(const IBackground &bg)
Definition: Simulation.cpp:257
const Instrument & instrument() const
Definition: Simulation.h:55
std::unique_ptr< IBackground > m_background
Definition: Simulation.h:157
SampleProvider m_sample_provider
Definition: Simulation.h:154
DistributionHandler m_distribution_handler
Definition: Simulation.h:155
ProgressHandler m_progress
Definition: Simulation.h:153

References initialize(), m_background, and setBackground().

Here is the call graph for this function:

Member Function Documentation

◆ clone()

virtual Simulation* Simulation::clone ( ) const
pure virtual

◆ prepareSimulation()

void Simulation::prepareSimulation ( )
virtual

Put into a clean state for running a simulation.

Reimplemented in SpecularSimulation, Simulation2D, GISASSimulation, and OffSpecSimulation.

Definition at line 189 of file Simulation.cpp.

190 {
193  throw std::runtime_error(
194  "Error in Simulation::prepareSimulation(): non-default materials of"
195  " several different types are used in the sample provided");
196  gsl_set_error_handler_off();
197 }
void updateSample()
Generates new sample if sample builder defined.
const MultiLayer * sample() const
Returns current sample.
bool ContainsCompatibleMaterials(const MultiLayer &multilayer)
Returns true if the multilayer contains non-default materials of one type only.

References MultiLayerUtils::ContainsCompatibleMaterials(), m_sample_provider, SampleProvider::sample(), and SampleProvider::updateSample().

Referenced by Simulation2D::prepareSimulation(), SpecularSimulation::prepareSimulation(), and runSimulation().

Here is the call graph for this function:

◆ runSimulation()

void Simulation::runSimulation ( )

Run a simulation, possibly averaged over parameter distributions.

Run simulation with possible averaging over parameter distributions.

Definition at line 200 of file Simulation.cpp.

201 {
203 
204  const size_t total_size = numberOfSimulationElements();
205  size_t param_combinations = m_distribution_handler.getTotalNumberOfSamples();
206 
207  m_progress.reset();
208  m_progress.setExpectedNTicks(param_combinations * total_size);
209 
210  // restrict calculation to current batch
211  const size_t n_batches = m_options.getNumberOfBatches();
212  const size_t current_batch = m_options.getCurrentBatch();
213 
214  const size_t batch_start = getStartIndex(n_batches, current_batch, total_size);
215  const size_t batch_size = getNumberOfElements(n_batches, current_batch, total_size);
216  if (batch_size == 0)
217  return;
218 
219  std::unique_ptr<ParameterPool> param_pool(createParameterTree());
220  for (size_t index = 0; index < param_combinations; ++index) {
221  double weight = m_distribution_handler.setParameterValues(param_pool.get(), index);
222  runSingleSimulation(batch_start, batch_size, weight);
223  }
227 }
double setParameterValues(ParameterPool *p_parameter_pool, size_t index)
set the parameter values of the simulation object to a specific combination of values,...
void setParameterToMeans(ParameterPool *p_parameter_pool) const
Sets mean distribution values to the parameter pool.
size_t getTotalNumberOfSamples() const
get the total number of parameter value combinations (product of the individual sizes of each paramet...
ParameterPool * createParameterTree() const
Creates new parameter pool, with all local parameters and those of its children.
Definition: INode.cpp:116
void setExpectedNTicks(size_t n)
unsigned getNumberOfBatches() const
unsigned getCurrentBatch() const
virtual void moveDataFromCache()=0
virtual void prepareSimulation()
Put into a clean state for running a simulation.
Definition: Simulation.cpp:189
virtual void transferResultsToIntensityMap()
Creates the appropriate data structure (e.g.
Definition: Simulation.h:110
void runSingleSimulation(size_t batch_start, size_t batch_size, double weight=1.0)
Runs a single simulation with fixed parameter values.
Definition: Simulation.cpp:289
virtual size_t numberOfSimulationElements() const =0
Gets the number of elements this simulation needs to calculate.
size_t getStartIndex(size_t n_handlers, size_t current_handler, size_t n_elements)
Definition: Simulation.cpp:53
size_t getNumberOfElements(size_t n_handlers, size_t current_handler, size_t n_elements)
Definition: Simulation.cpp:62

References INode::createParameterTree(), SimulationOptions::getCurrentBatch(), SimulationOptions::getNumberOfBatches(), anonymous_namespace{Simulation.cpp}::getNumberOfElements(), anonymous_namespace{Simulation.cpp}::getStartIndex(), DistributionHandler::getTotalNumberOfSamples(), m_distribution_handler, m_options, m_progress, moveDataFromCache(), numberOfSimulationElements(), prepareSimulation(), ProgressHandler::reset(), runSingleSimulation(), ProgressHandler::setExpectedNTicks(), DistributionHandler::setParameterToMeans(), DistributionHandler::setParameterValues(), and transferResultsToIntensityMap().

Here is the call graph for this function:

◆ runMPISimulation()

void Simulation::runMPISimulation ( )

Run a simulation in a MPI environment.

Definition at line 229 of file Simulation.cpp.

230 {
231  MPISimulation ompi;
232  ompi.runSimulation(this);
233 }
void runSimulation(Simulation *simulation)

References MPISimulation::runSimulation().

Here is the call graph for this function:

◆ setInstrument()

void Simulation::setInstrument ( const Instrument instrument_)

Definition at line 235 of file Simulation.cpp.

236 {
237  m_instrument = instrument_;
239 }
virtual void updateIntensityMap()
Definition: Simulation.h:115

References m_instrument, and updateIntensityMap().

Here is the call graph for this function:

◆ instrument() [1/2]

const Instrument& Simulation::instrument ( ) const
inline

Definition at line 55 of file Simulation.h.

55 { return m_instrument; }

References m_instrument.

Referenced by Simulation2D::addMask(), OffSpecSimulation::checkInitialization(), convertData(), DepthProbeSimulation::createUnitConverter(), OffSpecSimulation::createUnitConverter(), SimulationToPython::defineDetector(), SimulationToPython::defineDetectorPolarizationAnalysis(), SimulationToPython::defineDetectorResolutionFunction(), SimulationToPython::defineGISASBeam(), SimulationToPython::defineMasks(), SimulationToPython::defineOffSpecBeam(), SimulationToPython::defineSpecularScan(), Simulation2D::generateSimulationElements(), getBeamIntensity(), getChildren(), DepthProbeSimulation::initialize(), SpecularSimulation::initialize(), DepthProbeSimulation::initSimulationElementVector(), GISASSimulation::initSimulationElementVector(), OffSpecSimulation::initSimulationElementVector(), SpecularSimulation::initSimulationElementVector(), GISASSimulation::intensityMapSize(), OffSpecSimulation::intensityMapSize(), Simulation2D::maskAll(), DepthProbeSimulation::normalize(), SpecularSimulation::normalize(), GISASSimulation::prepareSimulation(), Simulation2D::prepareSimulation(), SpecularSimulation::prepareSimulation(), removeDetectorResolutionFunction(), Simulation2D::removeMasks(), GISASSimulation::result(), OffSpecSimulation::result(), setAnalyzerProperties(), setBeamIntensity(), DepthProbeSimulation::setBeamParameters(), OffSpecSimulation::setBeamParameters(), GISASSimulation::setBeamParameters(), setBeamPolarization(), Simulation2D::setDetector(), Simulation2D::setDetectorParameters(), setDetectorResolutionFunction(), Simulation2D::setRegionOfInterest(), SpecularSimulation::setScan(), OffSpecSimulation::transferDetectorImage(), OffSpecSimulation::transferResultsToIntensityMap(), and OffSpecSimulation::updateIntensityMap().

◆ instrument() [2/2]

Instrument& Simulation::instrument ( )
inline

Definition at line 56 of file Simulation.h.

56 { return m_instrument; }

References m_instrument.

◆ setBeamIntensity()

void Simulation::setBeamIntensity ( double  intensity)

Definition at line 173 of file Simulation.cpp.

174 {
175  instrument().setBeamIntensity(intensity);
176 }
void setBeamIntensity(double intensity)
Definition: Instrument.cpp:106

References instrument(), and Instrument::setBeamIntensity().

Referenced by StandardSimulations::GISASWithMasks().

Here is the call graph for this function:

◆ getBeamIntensity()

double Simulation::getBeamIntensity ( ) const

Definition at line 178 of file Simulation.cpp.

179 {
180  return instrument().getBeamIntensity();
181 }
double getBeamIntensity() const
Definition: Instrument.cpp:116

References Instrument::getBeamIntensity(), and instrument().

Referenced by DepthProbeSimulation::normalize(), Simulation2D::normalize(), and SpecularSimulation::normalize().

Here is the call graph for this function:

◆ setBeamPolarization()

void Simulation::setBeamPolarization ( const kvector_t  bloch_vector)

Sets the beam polarization according to the given Bloch vector.

Definition at line 184 of file Simulation.cpp.

185 {
186  instrument().setBeamPolarization(bloch_vector);
187 }
void setBeamPolarization(const kvector_t bloch_vector)
Sets the beam's polarization according to the given Bloch vector.
Definition: Instrument.cpp:111

References instrument(), and Instrument::setBeamPolarization().

Referenced by StandardSimulations::BasicGISAS00(), StandardSimulations::BasicPolarizedGISAS(), StandardSimulations::MaxiGISAS00(), StandardSimulations::MiniGISASPolarizationMM(), StandardSimulations::MiniGISASPolarizationMP(), StandardSimulations::MiniGISASPolarizationPM(), and StandardSimulations::MiniGISASPolarizationPP().

Here is the call graph for this function:

◆ setDetectorResolutionFunction()

void Simulation::setDetectorResolutionFunction ( const IResolutionFunction2D resolution_function)

Definition at line 156 of file Simulation.cpp.

157 {
158  instrument().setDetectorResolutionFunction(resolution_function);
159 }
void setDetectorResolutionFunction(const IResolutionFunction2D &p_resolution_function)
Sets detector resolution function.
Definition: Instrument.cpp:72

References instrument(), and Instrument::setDetectorResolutionFunction().

Referenced by StandardSimulations::MiniGISASDetectorResolution().

Here is the call graph for this function:

◆ removeDetectorResolutionFunction()

void Simulation::removeDetectorResolutionFunction ( )

Definition at line 161 of file Simulation.cpp.

162 {
164 }
void removeDetectorResolution()
Removes detector resolution function.
Definition: Instrument.cpp:77

References instrument(), and Instrument::removeDetectorResolution().

Here is the call graph for this function:

◆ setAnalyzerProperties()

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

Sets the polarization analyzer characteristics of the detector.

Definition at line 167 of file Simulation.cpp.

169 {
170  instrument().setAnalyzerProperties(direction, efficiency, total_transmission);
171 }
void setAnalyzerProperties(const kvector_t direction, double efficiency, double total_transmission)
Sets the polarization analyzer characteristics of the detector.
Definition: Instrument.cpp:167

References instrument(), and Instrument::setAnalyzerProperties().

Referenced by StandardSimulations::BasicGISAS00(), StandardSimulations::BasicPolarizedGISAS(), StandardSimulations::MaxiGISAS00(), StandardSimulations::MiniGISASPolarizationMM(), StandardSimulations::MiniGISASPolarizationMP(), StandardSimulations::MiniGISASPolarizationPM(), and StandardSimulations::MiniGISASPolarizationPP().

Here is the call graph for this function:

◆ setSample()

void Simulation::setSample ( const MultiLayer sample)

The MultiLayer object will not be owned by the Simulation object.

Definition at line 242 of file Simulation.cpp.

243 {
245 }
void setSample(const MultiLayer &multilayer)
const MultiLayer * sample() const
Definition: Simulation.cpp:247

References m_sample_provider, sample(), and SampleProvider::setSample().

Here is the call graph for this function:

◆ sample()

const MultiLayer * Simulation::sample ( ) const

Definition at line 247 of file Simulation.cpp.

248 {
249  return m_sample_provider.sample();
250 }

References m_sample_provider, and SampleProvider::sample().

Referenced by SimulationToPython::generateSimulationCode(), DepthProbeSimulation::generateSingleThreadedComputation(), Simulation2D::generateSingleThreadedComputation(), SpecularSimulation::generateSingleThreadedComputation(), setSample(), and DepthProbeSimulation::validityCheck().

Here is the call graph for this function:

◆ setSampleBuilder()

void Simulation::setSampleBuilder ( const std::shared_ptr< ISampleBuilder > &  sample_builder)

Definition at line 252 of file Simulation.cpp.

253 {
254  m_sample_provider.setBuilder(sample_builder);
255 }
void setBuilder(const std::shared_ptr< ISampleBuilder > &sample_builder)

References m_sample_provider, and SampleProvider::setBuilder().

Here is the call graph for this function:

◆ setBackground()

void Simulation::setBackground ( const IBackground bg)

Definition at line 257 of file Simulation.cpp.

258 {
259  m_background.reset(bg.clone());
261 }
virtual IBackground * clone() const =0
void registerChild(INode *node)
Definition: INode.cpp:58

References IBackground::clone(), m_background, and INode::registerChild().

Referenced by StandardSimulations::ConstantBackgroundGISAS(), and Simulation().

Here is the call graph for this function:

◆ background()

◆ intensityMapSize()

virtual size_t Simulation::intensityMapSize ( ) const
pure virtual

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

Implemented in SpecularSimulation, OffSpecSimulation, GISASSimulation, and DepthProbeSimulation.

◆ result()

virtual SimulationResult Simulation::result ( ) const
pure virtual

Returns the results of the simulation in a format that supports unit conversion and export to numpy arrays.

Implemented in SpecularSimulation, OffSpecSimulation, GISASSimulation, and DepthProbeSimulation.

Referenced by Simulation2D::generateSimulationElements(), getChildren(), and Simulation2D::rawResults().

◆ addParameterDistribution() [1/2]

void Simulation::addParameterDistribution ( const std::string &  param_name,
const IDistribution1D distribution,
size_t  nbr_samples,
double  sigma_factor = 0.0,
const RealLimits limits = RealLimits() 
)

Definition at line 273 of file Simulation.cpp.

276 {
277  ParameterDistribution par_distr(param_name, distribution, nbr_samples, sigma_factor, limits);
278  addParameterDistribution(par_distr);
279 }
A parametric distribution function, for use with any model parameter.
void addParameterDistribution(const std::string &param_name, const IDistribution1D &distribution, size_t nbr_samples, double sigma_factor=0.0, const RealLimits &limits=RealLimits())
Definition: Simulation.cpp:273

Referenced by StandardSimulations::MiniGISASBeamDivergence().

◆ addParameterDistribution() [2/2]

void Simulation::addParameterDistribution ( const ParameterDistribution par_distr)

Definition at line 281 of file Simulation.cpp.

282 {
283  validateParametrization(par_distr);
285 }
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: Simulation.h:135

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

Here is the call graph for this function:

◆ getDistributionHandler()

const DistributionHandler& Simulation::getDistributionHandler ( ) const
inline

Definition at line 89 of file Simulation.h.

89 { return m_distribution_handler; }

References m_distribution_handler.

Referenced by SimulationToPython::defineParameterDistributions().

◆ setOptions()

void Simulation::setOptions ( const SimulationOptions options)
inline

Definition at line 91 of file Simulation.h.

91 { m_options = options; }
const SimulationOptions & options() const
Definition: Simulation.h:120

References m_options, and options().

Here is the call graph for this function:

◆ getOptions() [1/2]

const SimulationOptions& Simulation::getOptions ( ) const
inline

◆ getOptions() [2/2]

SimulationOptions& Simulation::getOptions ( )
inline

Definition at line 93 of file Simulation.h.

93 { return m_options; }

References m_options.

◆ subscribe()

void Simulation::subscribe ( ProgressHandler::Callback_t  inform)
inline

Definition at line 95 of file Simulation.h.

95 { m_progress.subscribe(inform); }
void subscribe(ProgressHandler::Callback_t callback)

References m_progress, and ProgressHandler::subscribe().

Here is the call graph for this function:

◆ setTerminalProgressMonitor()

void Simulation::setTerminalProgressMonitor ( )

Initializes a progress monitor that prints to stdout.

Definition at line 145 of file Simulation.cpp.

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

References m_progress, and ProgressHandler::subscribe().

Here is the call graph for this function:

◆ getChildren()

std::vector< const INode * > Simulation::getChildren ( ) const
virtual

Returns a vector of children (const).

Reimplemented from INode.

Definition at line 263 of file Simulation.cpp.

264 {
265  std::vector<const INode*> result;
266  result.push_back(&instrument());
268  if (m_background)
269  result.push_back(m_background.get());
270  return result;
271 }
std::vector< const INode * > getChildren() const override
Returns a vector of children (const).
virtual SimulationResult result() const =0
Returns the results of the simulation in a format that supports unit conversion and export to numpy a...

References SampleProvider::getChildren(), instrument(), m_background, m_sample_provider, and result().

Here is the call graph for this function:

◆ convertData()

SimulationResult Simulation::convertData ( const OutputData< double > &  data,
bool  put_masked_areas_to_zero = true 
)

Convert user data to SimulationResult object for later drawing in various axes units.

User data will be cropped to the ROI defined in the simulation, amplitudes in areas corresponding to the masked areas of the detector will be set to zero.

Definition at line 317 of file Simulation.cpp.

319 {
320  auto converter = UnitConverterUtils::createConverter(*this);
321  auto roi_data = UnitConverterUtils::createOutputData(*converter, converter->defaultUnits());
322 
323  const IDetector& detector = instrument().detector();
324 
325  if (roi_data->hasSameDimensions(data)) {
326  // data is already cropped to ROI
327  if (put_masked_areas_to_zero) {
328  detector.iterate(
329  [&](IDetector::const_iterator it) {
330  (*roi_data)[it.roiIndex()] = data[it.roiIndex()];
331  },
332  /*visit_masked*/ false);
333  } else {
334  roi_data->setRawDataVector(data.getRawDataVector());
335  }
336 
337  } else if (detHasSameDimensions(detector, data)) {
338  // exp data has same shape as the detector, we have to put orig data to smaller roi map
339  detector.iterate(
340  [&](IDetector::const_iterator it) {
341  (*roi_data)[it.roiIndex()] = data[it.detectorIndex()];
342  },
343  /*visit_masked*/ !put_masked_areas_to_zero);
344 
345  } else {
346  throw std::runtime_error("FitObject::init_dataset() -> Error. Detector and exp data have "
347  "different shape.");
348  }
349 
350  return SimulationResult(*roi_data, *converter);
351 }
Abstract detector interface.
Definition: IDetector.h:36
void iterate(std::function< void(const_iterator)> func, bool visit_masks=false) const
Definition: IDetector.cpp:201
IDetector & detector()
Definition: Instrument.cpp:133
std::vector< T > getRawDataVector() const
Returns copy of raw data vector.
Definition: OutputData.h:335
An iterator for SimulationArea.
Wrapper around OutputData<double> that also provides unit conversions.
std::unique_ptr< IUnitConverter > createConverter(const Simulation &simulation)
std::unique_ptr< OutputData< double > > createOutputData(const IUnitConverter &converter, Axes::Units units)
Returns zero-valued output data array in specified units.
bool detHasSameDimensions(const IDetector &detector, const OutputData< double > &data)
Definition: Simulation.cpp:33

References UnitConverterUtils::createConverter(), UnitConverterUtils::createOutputData(), Instrument::detector(), SimulationAreaIterator::detectorIndex(), anonymous_namespace{Simulation.cpp}::detHasSameDimensions(), OutputData< T >::getRawDataVector(), instrument(), IDetector::iterate(), and SimulationAreaIterator::roiIndex().

Here is the call graph for this function:

◆ transferResultsToIntensityMap()

virtual void Simulation::transferResultsToIntensityMap ( )
inlineprotectedvirtual

Creates the appropriate data structure (e.g.

2D intensity map) from the calculated SimulationElement objects

Reimplemented in OffSpecSimulation.

Definition at line 110 of file Simulation.h.

110 {}

Referenced by runSimulation(), Simulation2D::setRawResults(), and SpecularSimulation::setRawResults().

◆ initSimulationElementVector()

virtual void Simulation::initSimulationElementVector ( )
protectedpure virtual

Initializes the vector of Simulation elements.

Implemented in SpecularSimulation, OffSpecSimulation, GISASSimulation, and DepthProbeSimulation.

Referenced by runSingleSimulation(), and Simulation2D::setRawResults().

◆ updateIntensityMap()

virtual void Simulation::updateIntensityMap ( )
inlineprotectedvirtual

Reimplemented in OffSpecSimulation.

Definition at line 115 of file Simulation.h.

115 {}

Referenced by Simulation2D::setDetectorParameters(), and setInstrument().

◆ numberOfSimulationElements()

virtual size_t Simulation::numberOfSimulationElements ( ) const
protectedpure virtual

Gets the number of elements this simulation needs to calculate.

Implemented in SpecularSimulation, Simulation2D, DepthProbeSimulation, and OffSpecSimulation.

Referenced by runSimulation().

◆ options()

◆ progress()

◆ initialize()

void Simulation::initialize ( )
private

Definition at line 138 of file Simulation.cpp.

139 {
142 }

References m_instrument, m_sample_provider, and INode::registerChild().

Referenced by Simulation().

Here is the call graph for this function:

◆ runSingleSimulation()

void Simulation::runSingleSimulation ( size_t  batch_start,
size_t  batch_size,
double  weight = 1.0 
)
private

Runs a single simulation with fixed parameter values.

If desired, the simulation is run in several threads.

Definition at line 289 of file Simulation.cpp.

290 {
292 
293  const size_t n_threads = m_options.getNumberOfThreads();
294  ASSERT(n_threads > 0);
295 
296  std::vector<std::unique_ptr<IComputation>> computations;
297 
298  for (size_t i_thread = 0; i_thread < n_threads;
299  ++i_thread) { // Distribute computations by threads
300  const size_t thread_start = batch_start + getStartIndex(n_threads, i_thread, batch_size);
301  const size_t thread_size = getNumberOfElements(n_threads, i_thread, batch_size);
302  if (thread_size == 0)
303  break;
304  computations.push_back(generateSingleThreadedComputation(thread_start, thread_size));
305  }
306  runComputations(std::move(computations));
307 
308  normalize(batch_start, batch_size);
309  addBackgroundIntensity(batch_start, batch_size);
310  addDataToCache(weight);
311 }
#define ASSERT(condition)
Definition: Assert.h:26
unsigned getNumberOfThreads() const
virtual std::unique_ptr< IComputation > generateSingleThreadedComputation(size_t start, size_t n_elements)=0
Generate a single threaded computation for a given range of simulation elements.
virtual void initSimulationElementVector()=0
Initializes the vector of Simulation elements.
virtual void addDataToCache(double weight)=0
virtual void addBackgroundIntensity(size_t start_ind, size_t n_elements)=0
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.
void runComputations(std::vector< std::unique_ptr< IComputation >> computations)
Definition: Simulation.cpp:71

References addBackgroundIntensity(), addDataToCache(), ASSERT, generateSingleThreadedComputation(), anonymous_namespace{Simulation.cpp}::getNumberOfElements(), SimulationOptions::getNumberOfThreads(), anonymous_namespace{Simulation.cpp}::getStartIndex(), initSimulationElementVector(), m_options, normalize(), and anonymous_namespace{Simulation.cpp}::runComputations().

Referenced by runSimulation().

Here is the call graph for this function:

◆ generateSingleThreadedComputation()

virtual std::unique_ptr<IComputation> Simulation::generateSingleThreadedComputation ( size_t  start,
size_t  n_elements 
)
privatepure virtual

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

Parameters
startIndex of the first element to include into computation
n_elementsNumber of elements to process

Implemented in SpecularSimulation, Simulation2D, and DepthProbeSimulation.

Referenced by runSingleSimulation().

◆ validateParametrization()

virtual void Simulation::validateParametrization ( const ParameterDistribution ) const
inlineprivatevirtual

Checks the distribution validity for simulation.

Reimplemented in SpecularSimulation, OffSpecSimulation, and DepthProbeSimulation.

Definition at line 135 of file Simulation.h.

135 {}

Referenced by addParameterDistribution().

◆ addBackgroundIntensity()

virtual void Simulation::addBackgroundIntensity ( size_t  start_ind,
size_t  n_elements 
)
privatepure virtual

◆ normalize()

virtual void Simulation::normalize ( size_t  start_ind,
size_t  n_elements 
)
privatepure virtual

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

Implemented in SpecularSimulation, Simulation2D, and DepthProbeSimulation.

Referenced by runSingleSimulation().

◆ addDataToCache()

virtual void Simulation::addDataToCache ( double  weight)
privatepure virtual

◆ moveDataFromCache()

virtual void Simulation::moveDataFromCache ( )
privatepure virtual

◆ rawResults()

virtual std::vector<double> Simulation::rawResults ( ) const
privatepure virtual

◆ setRawResults()

virtual void Simulation::setRawResults ( const std::vector< double > &  raw_data)
privatepure virtual

◆ transferToCPP()

virtual void ICloneable::transferToCPP ( )
inlinevirtualinherited

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

Definition at line 34 of file ICloneable.h.

◆ accept()

virtual void INode::accept ( INodeVisitor visitor) const
pure virtualinherited

Calls the INodeVisitor's visit method.

Implemented in FormFactorSphereLogNormalRadius, FormFactorSphereGaussianRadius, FormFactorGaussSphere, FormFactorDecoratorRotation, FormFactorDecoratorPositionFactor, FormFactorDecoratorMaterial, ParticleDistribution, ParticleCoreShell, ParticleComposition, Particle, MesoCrystal, FormFactorWeighted, FormFactorCrystal, FormFactorCoreShell, Crystal, Layer, FormFactorTruncatedSpheroid, FormFactorTruncatedSphere, FormFactorTruncatedCube, FormFactorTetrahedron, FormFactorSawtoothRippleLorentz, FormFactorSawtoothRippleGauss, FormFactorSawtoothRippleBox, FormFactorPyramid, FormFactorPrism6, FormFactorPrism3, FormFactorLongBoxLorentz, FormFactorLongBoxGauss, FormFactorIcosahedron, FormFactorHollowSphere, FormFactorHemiEllipsoid, FormFactorFullSpheroid, FormFactorFullSphere, FormFactorEllipsoidalCylinder, FormFactorDot, FormFactorDodecahedron, FormFactorCylinder, FormFactorCuboctahedron, FormFactorCosineRippleLorentz, FormFactorCosineRippleGauss, FormFactorCosineRippleBox, FormFactorCone6, FormFactorCone, FormFactorCantellatedCube, FormFactorBox, FormFactorAnisoPyramid, FTDistribution1DVoigt, FTDistribution1DCosine, FTDistribution1DTriangle, FTDistribution1DGate, FTDistribution1DGauss, FTDistribution1DCauchy, InterferenceFunctionTwin, InterferenceFunctionRadialParaCrystal, InterferenceFunctionNone, InterferenceFunctionHardDisk, InterferenceFunctionFinite3DLattice, InterferenceFunctionFinite2DLattice, InterferenceFunction3DLattice, InterferenceFunction2DSuperLattice, InterferenceFunction2DParaCrystal, InterferenceFunction2DLattice, InterferenceFunction1DLattice, SpecularSimulation, DepthProbeSimulation, FormFactorDWBAPol, FormFactorDWBA, FormFactorBAPol, Lattice, MisesGaussPeakShape, MisesFisherGaussPeakShape, LorentzFisherPeakShape, GaussFisherPeakShape, IsotropicLorentzPeakShape, IsotropicGaussPeakShape, SphericalDetector, SpecularDetector1D, FootprintSquare, FootprintGauss, Beam, GISASSimulation, PoissonNoiseBackground, ConstantBackground, MultiLayer, ParticleLayout, SampleProvider, SampleBuilderNode, HexagonalLattice, SquareLattice, BasicLattice, FormFactorBarLorentz, FormFactorBarGauss, FTDistribution2DVoigt, FTDistribution2DCone, FTDistribution2DGate, FTDistribution2DGauss, FTDistribution2DCauchy, FTDecayFunction2DVoigt, FTDecayFunction2DGauss, FTDecayFunction2DCauchy, FTDecayFunction1DVoigt, FTDecayFunction1DTriangle, FTDecayFunction1DGauss, FTDecayFunction1DCauchy, DistributionTrapezoid, DistributionCosine, DistributionLogNormal, DistributionGaussian, DistributionLorentz, DistributionGate, ResolutionFunction2DGaussian, ConvolutionDetectorResolution, Instrument, RectangularDetector, IsGISAXSDetector, DetectionProperties, OffSpecSimulation, ILayout, LayerRoughness, LayerInterface, RotationEuler, RotationZ, RotationY, RotationX, IdentityRotation, and IAbstractParticle.

Referenced by VisitNodesPostorder(), and VisitNodesPreorder().

◆ treeToString()

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

Returns multiline string representing tree structure below the node.

Definition at line 53 of file INode.cpp.

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

References NodeUtils::nodeToString().

Here is the call graph for this function:

◆ registerChild()

void INode::registerChild ( INode node)
inherited

Definition at line 58 of file INode.cpp.

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

References ASSERT, and INode::setParent().

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

Here is the call graph for this function:

◆ setParent()

void INode::setParent ( const INode newParent)
virtualinherited

Reimplemented in SampleProvider.

Definition at line 69 of file INode.cpp.

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

References INode::m_parent.

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

◆ parent() [1/2]

const INode * INode::parent ( ) const
inherited

◆ parent() [2/2]

INode * INode::parent ( )
inherited

Definition at line 79 of file INode.cpp.

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

References INode::m_parent.

◆ copyNumber()

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

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

Definition at line 84 of file INode.cpp.

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

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

Referenced by INode::displayName().

Here is the call graph for this function:

◆ displayName()

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

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

Definition at line 105 of file INode.cpp.

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

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

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

Here is the call graph for this function:

◆ createParameterTree()

ParameterPool * INode::createParameterTree ( ) const
virtualinherited

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

Reimplemented from IParameterized.

Definition at line 116 of file INode.cpp.

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

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

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

Here is the call graph for this function:

◆ parameterPool()

ParameterPool* IParameterized::parameterPool ( ) const
inlineinherited

Returns pointer to the parameter pool.

Definition at line 38 of file IParameterized.h.

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

References IParameterized::m_pool.

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

◆ parametersToString()

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

Returns multiline string representing available parameters.

Definition at line 40 of file IParameterized.cpp.

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

References IParameterized::createParameterTree().

Here is the call graph for this function:

◆ registerParameter()

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

Definition at line 48 of file IParameterized.cpp.

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

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

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

Here is the call graph for this function:

◆ registerVector()

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

Definition at line 54 of file IParameterized.cpp.

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

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

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

Here is the call graph for this function:

◆ setParameterValue()

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

Definition at line 62 of file IParameterized.cpp.

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

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

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

Here is the call graph for this function:

◆ setVectorValue()

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

Definition at line 75 of file IParameterized.cpp.

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

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

Here is the call graph for this function:

◆ parameter()

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

◆ onChange()

◆ removeParameter()

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

◆ removeVector()

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

Definition at line 93 of file IParameterized.cpp.

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

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

Referenced by IParticle::registerPosition().

Here is the call graph for this function:

◆ XComponentName()

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

◆ YComponentName()

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

Definition at line 105 of file IParameterized.cpp.

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

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

◆ ZComponentName()

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

Definition at line 110 of file IParameterized.cpp.

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

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

◆ setName()

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

Definition at line 68 of file IParameterized.h.

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

References IParameterized::m_name.

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

◆ getName()

Friends And Related Function Documentation

◆ MPISimulation

friend class MPISimulation
friend

Definition at line 103 of file Simulation.h.

Member Data Documentation

◆ m_options

SimulationOptions Simulation::m_options
private

Definition at line 152 of file Simulation.h.

Referenced by getOptions(), options(), runSimulation(), runSingleSimulation(), and setOptions().

◆ m_progress

ProgressHandler Simulation::m_progress
private

Definition at line 153 of file Simulation.h.

Referenced by progress(), runSimulation(), setTerminalProgressMonitor(), and subscribe().

◆ m_sample_provider

SampleProvider Simulation::m_sample_provider
private

◆ m_distribution_handler

DistributionHandler Simulation::m_distribution_handler
private

Definition at line 155 of file Simulation.h.

Referenced by addParameterDistribution(), getDistributionHandler(), and runSimulation().

◆ m_instrument

Instrument Simulation::m_instrument
private

Definition at line 156 of file Simulation.h.

Referenced by initialize(), instrument(), and setInstrument().

◆ m_background

std::unique_ptr<IBackground> Simulation::m_background
private

Definition at line 157 of file Simulation.h.

Referenced by background(), getChildren(), setBackground(), and Simulation().

◆ m_parent

const INode* INode::m_parent {nullptr}
privateinherited

Definition at line 81 of file INode.h.

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

◆ m_NP

const size_t INode::m_NP
protectedinherited

Definition at line 86 of file INode.h.

Referenced by INode::INode().

◆ m_P

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

Definition at line 87 of file INode.h.

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

◆ m_name

std::string IParameterized::m_name
privateinherited

Definition at line 72 of file IParameterized.h.

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

◆ m_pool

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

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