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

Public Member Functions

 FormFactorSphereGaussianRadius (const std::vector< double > P)
 
 FormFactorSphereGaussianRadius (double mean, double sigma)
 
FormFactorSphereGaussianRadiusclone () const override final
 
void accept (INodeVisitor *visitor) const override final
 
double radialExtension () const override final
 
complex_t evaluate_for_q (cvector_t q) const override final
 
void setAmbientMaterial (const Material &) override
 
complex_t evaluate (const WavevectorInfo &wavevectors) const override
 
Eigen::Matrix2cd evaluatePol (const WavevectorInfo &wavevectors) const override
 
virtual double bottomZ (const IRotation &rotation) const override
 
virtual double topZ (const IRotation &rotation) const override
 
IFormFactorcreateSlicedFormFactor (ZLimits limits, const IRotation &rot, kvector_t translation) const
 
virtual double volume () const
 
virtual void setSpecularInfo (std::unique_ptr< const ILayerRTCoefficients >, std::unique_ptr< const ILayerRTCoefficients >)
 
virtual const Materialmaterial () const
 
std::vector< const Material * > containedMaterials () const
 
virtual void transferToCPP ()
 
virtual std::string treeToString () const
 
void registerChild (INode *node)
 
virtual std::vector< const INode * > getChildren () const
 
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
 
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

void onChange () override final
 
bool canSliceAnalytically (const IRotation &rot) const override
 
virtual Eigen::Matrix2cd evaluate_for_q_pol (cvector_t q) const
 
SlicingEffects computeSlicingEffects (ZLimits limits, const kvector_t &position, double height) const
 
virtual IFormFactorsliceFormFactor (ZLimits limits, const IRotation &rot, kvector_t translation) const
 

Static Protected Member Functions

static double BottomZ (const std::vector< kvector_t > &vertices, const IRotation &rotation)
 
static double TopZ (const std::vector< kvector_t > &vertices, const IRotation &rotation)
 

Protected Attributes

std::unique_ptr< IShapemP_shape
 
const size_t m_NP
 
std::vector< double > m_P
 

Private Member Functions

double calculateMeanR3 () const
 

Private Attributes

const double & m_mean
 
const double & m_sigma
 
double m_mean_r3
 
const INodem_parent {nullptr}
 
std::string m_name
 
std::unique_ptr< ParameterPoolm_pool
 

Detailed Description

A sphere with gaussian radius distribution.

Definition at line 24 of file FormFactorSphereGaussianRadius.h.

Constructor & Destructor Documentation

◆ FormFactorSphereGaussianRadius() [1/2]

FormFactorSphereGaussianRadius::FormFactorSphereGaussianRadius ( const std::vector< double >  P)

Definition at line 20 of file FormFactorSphereGaussianRadius.cpp.

21  : IFormFactorBorn({"FormFactorSphereGaussianRadius",
22  "class_tooltip",
23  {{"MeanRadius", "nm", "para_tooltip", 0, +INF, 0},
24  {"SigmaRadius", "nm", "para_tooltip", 0, +INF, 0}}},
25  P),
26  m_mean(m_P[0]), m_sigma(m_P[1])
27 {
29  onChange();
30 }
const double INF
Definition: INode.h:24
const double & m_mean
This is the mean radius.
double m_mean_r3
This is the radius that gives the mean volume.
void onChange() override final
Action to be taken in inherited class when a parameter has changed.
std::vector< double > m_P
Definition: INode.h:87

References INF.

Referenced by clone().

◆ FormFactorSphereGaussianRadius() [2/2]

FormFactorSphereGaussianRadius::FormFactorSphereGaussianRadius ( double  mean,
double  sigma 
)

Definition at line 32 of file FormFactorSphereGaussianRadius.cpp.

33  : FormFactorSphereGaussianRadius(std::vector<double>{mean, sigma})
34 {
35 }
FormFactorSphereGaussianRadius(const std::vector< double > P)

Member Function Documentation

◆ clone()

FormFactorSphereGaussianRadius* FormFactorSphereGaussianRadius::clone ( ) const
inlinefinaloverridevirtual

Returns a clone of this ISample object.

Implements IFormFactorBorn.

Definition at line 30 of file FormFactorSphereGaussianRadius.h.

31  {
33  }

References FormFactorSphereGaussianRadius(), m_mean, and m_sigma.

Here is the call graph for this function:

◆ accept()

void FormFactorSphereGaussianRadius::accept ( INodeVisitor visitor) const
inlinefinaloverridevirtual

Calls the INodeVisitor's visit method.

Implements INode.

Definition at line 35 of file FormFactorSphereGaussianRadius.h.

35 { visitor->visit(this); }
virtual void visit(const BasicLattice *)
Definition: INodeVisitor.h:154

◆ radialExtension()

double FormFactorSphereGaussianRadius::radialExtension ( ) const
inlinefinaloverridevirtual

Returns the (approximate in some cases) radial size of the particle of this form factor's shape.

This is used for SSCA calculations

Implements IFormFactor.

Definition at line 37 of file FormFactorSphereGaussianRadius.h.

37 { return m_mean; }

References m_mean.

◆ evaluate_for_q()

complex_t FormFactorSphereGaussianRadius::evaluate_for_q ( cvector_t  q) const
finaloverridevirtual

Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.

This method is public only for convenience of plotting form factors in Python.

Implements IFormFactorBorn.

Definition at line 37 of file FormFactorSphereGaussianRadius.cpp.

38 {
39  double q2 = std::norm(q.x()) + std::norm(q.y()) + std::norm(q.z());
40  double dw = std::exp(-q2 * m_sigma * m_sigma / 2.0);
41  return dw * exp_I(q.z() * m_mean_r3) * someff::ffSphere(q, m_mean_r3);
42  // TODO: don't center at bottom; revise mesocrystal1.py
43 }
complex_t exp_I(complex_t z)
Returns exp(I*z), where I is the imaginary unit.
Definition: Complex.h:30
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
complex_t ffSphere(cvector_t q, double R)
Returns the form factor of a sphere of radius R.

References exp_I(), someff::ffSphere(), m_mean_r3, m_sigma, BasicVector3D< T >::x(), BasicVector3D< T >::y(), and BasicVector3D< T >::z().

Here is the call graph for this function:

◆ onChange()

void FormFactorSphereGaussianRadius::onChange ( )
finaloverrideprotectedvirtual

Action to be taken in inherited class when a parameter has changed.

Reimplemented from IParameterized.

Definition at line 45 of file FormFactorSphereGaussianRadius.cpp.

46 {
47  mP_shape = std::make_unique<TruncatedEllipsoid>(m_mean, m_mean, m_mean, 2.0 * m_mean, 0.0);
48 }
std::unique_ptr< IShape > mP_shape
IShape object, used to retrieve vertices (which may be approximate in the case of round shapes).

References m_mean, and IFormFactorBorn::mP_shape.

◆ calculateMeanR3()

double FormFactorSphereGaussianRadius::calculateMeanR3 ( ) const
private

Definition at line 50 of file FormFactorSphereGaussianRadius.cpp.

51 {
52  return std::pow(m_mean * (m_mean * m_mean + 3.0 * m_sigma * m_sigma), 1.0 / 3.0);
53 }

References m_mean, and m_sigma.

◆ setAmbientMaterial()

void IFormFactorBorn::setAmbientMaterial ( const Material )
inlineoverridevirtualinherited

Passes the material in which this particle is embedded.

Implements IFormFactor.

Definition at line 40 of file IFormFactorBorn.h.

40 {}

◆ evaluate()

complex_t IFormFactorBorn::evaluate ( const WavevectorInfo wavevectors) const
overridevirtualinherited

Returns scattering amplitude for complex wavevectors ki, kf.

Implements IFormFactor.

Definition at line 31 of file IFormFactorBorn.cpp.

32 {
33  return evaluate_for_q(wavevectors.getQ());
34 }
virtual complex_t evaluate_for_q(cvector_t q) const =0
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
cvector_t getQ() const

References IFormFactorBorn::evaluate_for_q(), and WavevectorInfo::getQ().

Here is the call graph for this function:

◆ evaluatePol()

Eigen::Matrix2cd IFormFactorBorn::evaluatePol ( const WavevectorInfo wavevectors) const
overridevirtualinherited

Returns scattering amplitude for matrix interactions.

Reimplemented from IFormFactor.

Definition at line 36 of file IFormFactorBorn.cpp.

37 {
38  return evaluate_for_q_pol(wavevectors.getQ());
39 }
virtual Eigen::Matrix2cd evaluate_for_q_pol(cvector_t q) const
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case of matrix interactio...

References IFormFactorBorn::evaluate_for_q_pol(), and WavevectorInfo::getQ().

Here is the call graph for this function:

◆ bottomZ()

double IFormFactorBorn::bottomZ ( const IRotation rotation) const
overridevirtualinherited

Returns the z-coordinate of the lowest point in this shape after a given rotation.

Implements IFormFactor.

Reimplemented in IFormFactorPrism, IFormFactorPolyhedron, FormFactorFullSphere, and FormFactorDot.

Definition at line 41 of file IFormFactorBorn.cpp.

42 {
43  if (!mP_shape)
44  return 0;
45  return BottomZ(mP_shape->vertices(), rotation);
46 }
static double BottomZ(const std::vector< kvector_t > &vertices, const IRotation &rotation)
Calculates the z-coordinate of the lowest vertex after rotation.

References IFormFactorBorn::BottomZ(), and IFormFactorBorn::mP_shape.

Here is the call graph for this function:

◆ topZ()

double IFormFactorBorn::topZ ( const IRotation rotation) const
overridevirtualinherited

Returns the z-coordinate of the lowest point in this shape after a given rotation.

Implements IFormFactor.

Reimplemented in IFormFactorPrism, IFormFactorPolyhedron, FormFactorFullSphere, and FormFactorDot.

Definition at line 48 of file IFormFactorBorn.cpp.

49 {
50  if (!mP_shape)
51  return 0;
52  return TopZ(mP_shape->vertices(), rotation);
53 }
static double TopZ(const std::vector< kvector_t > &vertices, const IRotation &rotation)
Calculates the z-coordinate of the highest vertex after rotation.

References IFormFactorBorn::mP_shape, and IFormFactorBorn::TopZ().

Here is the call graph for this function:

◆ canSliceAnalytically()

bool IFormFactorBorn::canSliceAnalytically ( const IRotation rot) const
overrideprotectedvirtualinherited

Default implementation only allows rotations along z-axis.

Reimplemented from IFormFactor.

Definition at line 55 of file IFormFactorBorn.cpp.

56 {
57  if (rot.zInvariant())
58  return true;
59  return false;
60 }
bool zInvariant() const
Definition: Rotations.cpp:68

References IRotation::zInvariant().

Here is the call graph for this function:

◆ evaluate_for_q_pol()

Eigen::Matrix2cd IFormFactorBorn::evaluate_for_q_pol ( cvector_t  q) const
protectedvirtualinherited

Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case of matrix interactions.

Default implementation calls evaluate_for_q(q) and multiplies with the unit matrix.

Definition at line 62 of file IFormFactorBorn.cpp.

63 {
64  return evaluate_for_q(q) * Eigen::Matrix2cd::Identity();
65 }

References IFormFactorBorn::evaluate_for_q().

Referenced by IFormFactorBorn::evaluatePol().

Here is the call graph for this function:

◆ computeSlicingEffects()

SlicingEffects IFormFactorBorn::computeSlicingEffects ( ZLimits  limits,
const kvector_t position,
double  height 
) const
protectedinherited

Helper method for slicing.

Definition at line 67 of file IFormFactorBorn.cpp.

69 {
70  kvector_t new_position(position);
71  double z_bottom = position.z();
72  double z_top = position.z() + height;
73  OneSidedLimit lower_limit = limits.lowerLimit();
74  OneSidedLimit upper_limit = limits.upperLimit();
75  if (!upper_limit.m_limitless && !lower_limit.m_limitless
76  && lower_limit.m_value > upper_limit.m_value)
77  throw std::runtime_error(getName()
78  + "::sliceFormFactor error: "
79  "upperlimit < lowerlimit.");
80  double dz_top = upper_limit.m_limitless ? -1 : z_top - upper_limit.m_value;
81  double dz_bottom = lower_limit.m_limitless ? -1 : lower_limit.m_value - z_bottom;
82  if (dz_top < 0 && dz_bottom < 0)
83  throw std::runtime_error(getName()
84  + "::sliceFormFactor error: "
85  "shape didn't need to be sliced.");
86  if (dz_bottom > height)
87  throw std::runtime_error(getName()
88  + "::sliceFormFactor error: "
89  "interface outside shape.");
90  if (dz_top > height)
91  throw std::runtime_error(getName()
92  + "::sliceFormFactor error: "
93  "interface outside shape.");
94  if (dz_bottom < 0)
95  dz_bottom = 0;
96  if (dz_top < 0)
97  dz_top = 0;
98  if (dz_bottom > 0)
99  new_position.setZ(lower_limit.m_value);
100  return {new_position, dz_bottom, dz_top};
101 }
const std::string & getName() const
OneSidedLimit lowerLimit() const
Definition: ZLimits.cpp:39
OneSidedLimit upperLimit() const
Definition: ZLimits.cpp:44
Helper class that represents a onesided limit.
Definition: ZLimits.h:30
double m_value
Definition: ZLimits.h:32
bool m_limitless
Definition: ZLimits.h:31

References IParameterized::getName(), anonymous_namespace{BoxCompositionBuilder.cpp}::height, ZLimits::lowerLimit(), OneSidedLimit::m_limitless, OneSidedLimit::m_value, BasicVector3D< T >::setZ(), ZLimits::upperLimit(), and BasicVector3D< T >::z().

Referenced by FormFactorAnisoPyramid::sliceFormFactor(), FormFactorBox::sliceFormFactor(), FormFactorCone::sliceFormFactor(), FormFactorCone6::sliceFormFactor(), FormFactorCuboctahedron::sliceFormFactor(), FormFactorCylinder::sliceFormFactor(), FormFactorEllipsoidalCylinder::sliceFormFactor(), FormFactorFullSphere::sliceFormFactor(), FormFactorFullSpheroid::sliceFormFactor(), FormFactorLongBoxGauss::sliceFormFactor(), FormFactorLongBoxLorentz::sliceFormFactor(), FormFactorPrism3::sliceFormFactor(), FormFactorPrism6::sliceFormFactor(), FormFactorPyramid::sliceFormFactor(), FormFactorTetrahedron::sliceFormFactor(), FormFactorTruncatedSphere::sliceFormFactor(), and FormFactorTruncatedSpheroid::sliceFormFactor().

Here is the call graph for this function:

◆ BottomZ()

double IFormFactorBorn::BottomZ ( const std::vector< kvector_t > &  vertices,
const IRotation rotation 
)
staticprotectedinherited

Calculates the z-coordinate of the lowest vertex after rotation.

Definition at line 103 of file IFormFactorBorn.cpp.

104 {
105  ASSERT(vertices.size());
106  return algo::min_value(
107  vertices.begin(), vertices.end(),
108  [&](const kvector_t& vertex) -> double { return rotation.transformed(vertex).z(); });
109 }
#define ASSERT(condition)
Definition: Assert.h:26
double min_value(const Iterator &begin, const Iterator &end, const Evaluator &evaluate)
Returns the minimum value of function evaluate as applied to the elements of an iterator range.
Definition: Algorithms.h:54

References ASSERT, and algo::min_value().

Referenced by IFormFactorBorn::bottomZ(), IFormFactorPolyhedron::bottomZ(), and IFormFactorPrism::bottomZ().

Here is the call graph for this function:

◆ TopZ()

double IFormFactorBorn::TopZ ( const std::vector< kvector_t > &  vertices,
const IRotation rotation 
)
staticprotectedinherited

Calculates the z-coordinate of the highest vertex after rotation.

Definition at line 111 of file IFormFactorBorn.cpp.

112 {
113  ASSERT(vertices.size());
114  return algo::max_value(
115  vertices.begin(), vertices.end(),
116  [&](const kvector_t& vertex) -> double { return rotation.transformed(vertex).z(); });
117 }
double max_value(const Iterator &begin, const Iterator &end, const Evaluator &evaluate)
Returns the maximum value of function evaluate as applied to the elements of an iterator range.
Definition: Algorithms.h:65

References ASSERT, and algo::max_value().

Referenced by IFormFactorBorn::topZ(), IFormFactorPolyhedron::topZ(), and IFormFactorPrism::topZ().

Here is the call graph for this function:

◆ createSlicedFormFactor()

IFormFactor * IFormFactor::createSlicedFormFactor ( ZLimits  limits,
const IRotation rot,
kvector_t  translation 
) const
inherited

Creates a (possibly sliced) form factor with the given rotation and translation.

Definition at line 35 of file IFormFactor.cpp.

37 {
38  if (ShapeIsContainedInLimits(*this, limits, rot, translation))
39  return createTransformedFormFactor(*this, rot, translation);
40  if (ShapeOutsideLimits(*this, limits, rot, translation))
41  return nullptr;
42  if (canSliceAnalytically(rot))
43  return sliceFormFactor(limits, rot, translation);
44  throw std::runtime_error(getName()
45  + "::createSlicedFormFactor error: not supported for "
46  "the given rotation!");
47 }
IFormFactor * createTransformedFormFactor(const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:77
virtual IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const
Actually slices the form factor or throws an exception.
Definition: IFormFactor.cpp:72
virtual bool canSliceAnalytically(const IRotation &rot) const
Checks if slicing has a fast analytical solution.
Definition: IFormFactor.cpp:67
bool ShapeIsContainedInLimits(const IFormFactor &formfactor, ZLimits limits, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:94
bool ShapeOutsideLimits(const IFormFactor &formfactor, ZLimits limits, const IRotation &rot, kvector_t translation)

References IFormFactor::canSliceAnalytically(), createTransformedFormFactor(), IParameterized::getName(), anonymous_namespace{IFormFactor.cpp}::ShapeIsContainedInLimits(), anonymous_namespace{IFormFactor.cpp}::ShapeOutsideLimits(), and IFormFactor::sliceFormFactor().

Here is the call graph for this function:

◆ volume()

double IFormFactor::volume ( ) const
virtualinherited

Returns the total volume of the particle of this form factor's shape.

Reimplemented in FormFactorCrystal, IFormFactorPolyhedron, FormFactorBox, IFormFactorDecorator, FormFactorDWBAPol, FormFactorDWBA, FormFactorBAPol, and IFormFactorPrism.

Definition at line 56 of file IFormFactor.cpp.

57 {
58  auto zero_wavevectors = WavevectorInfo::GetZeroQ();
59  return std::abs(evaluate(zero_wavevectors));
60 }
virtual complex_t evaluate(const WavevectorInfo &wavevectors) const =0
Returns scattering amplitude for complex wavevectors ki, kf.
static WavevectorInfo GetZeroQ()

References IFormFactor::evaluate(), and WavevectorInfo::GetZeroQ().

Referenced by IFormFactorDecorator::volume(), and FormFactorCrystal::volume().

Here is the call graph for this function:

◆ setSpecularInfo()

void IFormFactor::setSpecularInfo ( std::unique_ptr< const ILayerRTCoefficients ,
std::unique_ptr< const ILayerRTCoefficients  
)
virtualinherited

Sets reflection/transmission info.

Reimplemented in FormFactorDWBAPol, and FormFactorDWBA.

Definition at line 62 of file IFormFactor.cpp.

64 {
65 }

◆ sliceFormFactor()

IFormFactor * IFormFactor::sliceFormFactor ( ZLimits  limits,
const IRotation rot,
kvector_t  translation 
) const
protectedvirtualinherited

Actually slices the form factor or throws an exception.

Reimplemented in FormFactorTruncatedSpheroid, FormFactorTruncatedSphere, FormFactorTetrahedron, FormFactorPyramid, FormFactorPrism6, FormFactorPrism3, FormFactorLongBoxLorentz, FormFactorLongBoxGauss, FormFactorFullSpheroid, FormFactorFullSphere, FormFactorEllipsoidalCylinder, FormFactorCylinder, FormFactorCuboctahedron, FormFactorCone6, FormFactorCone, FormFactorBox, and FormFactorAnisoPyramid.

Definition at line 72 of file IFormFactor.cpp.

73 {
74  throw std::runtime_error(getName() + "::sliceFormFactor error: not implemented!");
75 }

References IParameterized::getName().

Referenced by IFormFactor::createSlicedFormFactor().

Here is the call graph for this function:

◆ material()

◆ containedMaterials()

std::vector< const Material * > ISample::containedMaterials ( ) const
inherited

Returns set of unique materials contained in this ISample.

Definition at line 23 of file ISample.cpp.

24 {
25  std::vector<const Material*> result;
26  if (const Material* p_material = material())
27  result.push_back(p_material);
28  for (auto child : getChildren()) {
29  if (const ISample* sample = dynamic_cast<const ISample*>(child)) {
30  for (const Material* p_material : sample->containedMaterials())
31  result.push_back(p_material);
32  }
33  }
34  return result;
35 }
virtual std::vector< const INode * > getChildren() const
Returns a vector of children (const).
Definition: INode.cpp:64
Pure virtual base class for sample components and properties related to scattering.
Definition: ISample.h:28
virtual const Material * material() const
Returns nullptr, unless overwritten to return a specific material.
Definition: ISample.h:37
A wrapper for underlying material implementation.
Definition: Material.h:29

References INode::getChildren(), and ISample::material().

Referenced by MultiLayerUtils::ContainsCompatibleMaterials(), anonymous_namespace{ProcessedSample.cpp}::ContainsMagneticMaterial(), and SampleToPython::initLabels().

Here is the call graph for this function:

◆ transferToCPP()

virtual void ICloneable::transferToCPP ( )
inlinevirtualinherited

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

Definition at line 34 of file ICloneable.h.

◆ 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(), Simulation::initialize(), MesoCrystal::initialize(), Instrument::Instrument(), Beam::operator=(), Instrument::operator=(), Particle::Particle(), ParticleDistribution::ParticleDistribution(), IParticle::rotate(), ParticleLayout::setAndRegisterInterferenceFunction(), Simulation::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:

◆ getChildren()

◆ 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

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(), Simulation::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 }
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

◆ 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()

Member Data Documentation

◆ m_mean

const double& FormFactorSphereGaussianRadius::m_mean
private

This is the mean radius.

Definition at line 47 of file FormFactorSphereGaussianRadius.h.

Referenced by calculateMeanR3(), clone(), onChange(), and radialExtension().

◆ m_sigma

const double& FormFactorSphereGaussianRadius::m_sigma
private

Definition at line 48 of file FormFactorSphereGaussianRadius.h.

Referenced by calculateMeanR3(), clone(), and evaluate_for_q().

◆ m_mean_r3

double FormFactorSphereGaussianRadius::m_mean_r3
private

This is the radius that gives the mean volume.

Definition at line 49 of file FormFactorSphereGaussianRadius.h.

Referenced by evaluate_for_q().

◆ mP_shape

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