BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorDecoratorPositionFactor Class Reference

Decorates a form factor with a position dependent phase factor. More...

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

Public Member Functions

 FormFactorDecoratorPositionFactor (const IFormFactor &ff, const kvector_t &position)
 
void accept (INodeVisitor *visitor) const final
 Calls the INodeVisitor's visit method. More...
 
double bottomZ (const IRotation &rotation) const final
 Returns the z-coordinate of the lowest point in this shape after a given rotation. More...
 
FormFactorDecoratorPositionFactorclone () const final
 Returns a clone of this ISampleNode object. More...
 
std::vector< const Material * > containedMaterials () const
 Returns set of unique materials contained in this ISampleNode. More...
 
int copyNumber (const INode *node) const
 Returns copyNumber of child, which takes into account existence of children with same name. More...
 
ParameterPoolcreateParameterTree () const
 Creates new parameter pool, with all local parameters and those of its children. More...
 
IFormFactorcreateSlicedFormFactor (ZLimits limits, const IRotation &rot, kvector_t translation) const
 Creates a (possibly sliced) form factor with the given rotation and translation. More...
 
std::string displayName () const
 Returns display name, composed from the name of node and it's copy number. More...
 
complex_t evaluate (const WavevectorInfo &wavevectors) const final
 Returns scattering amplitude for complex wavevectors ki, kf. More...
 
Eigen::Matrix2cd evaluatePol (const WavevectorInfo &wavevectors) const final
 Returns scattering amplitude for matrix interactions. More...
 
virtual std::vector< const INode * > getChildren () const
 Returns a vector of children. More...
 
const IFormFactorgetFormFactor () const
 
const std::string & getName () const
 
bool isMagnetic () const
 Returns true if there is any magnetic material in this ISampleNode. More...
 
virtual const Materialmaterial () const
 Returns nullptr, unless overwritten to return a specific material. More...
 
virtual void onChange ()
 Action to be taken in inherited class when a parameter has changed. More...
 
RealParameterparameter (const std::string &name) const
 Returns parameter with given 'name'. More...
 
ParameterPoolparameterPool () const
 Returns pointer to the parameter pool. More...
 
std::string parametersToString () const
 Returns multiline string representing available parameters. More...
 
INodeparent ()
 
const INodeparent () const
 
std::vector< const INode * > progeny () const
 Returns a vector of all descendants. More...
 
double radialExtension () const override
 Returns the (approximate in some cases) radial size of the particle of this form factor's shape. More...
 
void registerChild (INode *node)
 
RealParameterregisterParameter (const std::string &name, double *parpointer)
 
void registerVector (const std::string &base_name, kvector_t *p_vec, const std::string &units="nm")
 
void removeParameter (const std::string &name)
 
void removeVector (const std::string &base_name)
 
void setAmbientMaterial (const Material &material) override
 Passes the material in which this particle is embedded. More...
 
void setName (const std::string &name)
 
void setParameterValue (const std::string &name, double value)
 
virtual void setParent (const INode *newParent)
 
virtual void setSpecularInfo (std::unique_ptr< const ILayerRTCoefficients >, std::unique_ptr< const ILayerRTCoefficients >)
 Sets reflection/transmission info. More...
 
void setVectorValue (const std::string &base_name, kvector_t value)
 
double topZ (const IRotation &rotation) const final
 Returns the z-coordinate of the lowest point in this shape after a given rotation. More...
 
virtual void transferToCPP ()
 Used for Python overriding of clone (see swig/tweaks.py) More...
 
virtual std::string treeToString () const
 Returns multiline string representing tree structure below the node. More...
 
double volume () const override
 Returns the total volume of the particle of this form factor's shape. More...
 

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

virtual bool canSliceAnalytically (const IRotation &rot) const
 Checks if slicing has a fast analytical solution. More...
 
virtual IFormFactorsliceFormFactor (ZLimits limits, const IRotation &rot, kvector_t translation) const
 Actually slices the form factor or throws an exception. More...
 

Static Protected Member Functions

static IFormFactorcreateTransformedFormFactor (const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
 

Protected Attributes

IFormFactorm_ff
 
const size_t m_NP
 
std::vector< double > m_P
 

Private Member Functions

complex_t getPositionFactor (const WavevectorInfo &wavevectors) const
 

Private Attributes

std::string m_name
 
const INodem_parent {nullptr}
 
std::unique_ptr< ParameterPoolm_pool
 parameter pool (kind of pointer-to-implementation) More...
 
kvector_t m_position
 

Detailed Description

Decorates a form factor with a position dependent phase factor.

Definition at line 28 of file FormFactorDecoratorPositionFactor.h.

Constructor & Destructor Documentation

◆ FormFactorDecoratorPositionFactor()

FormFactorDecoratorPositionFactor::FormFactorDecoratorPositionFactor ( const IFormFactor ff,
const kvector_t position 
)

Definition at line 19 of file FormFactorDecoratorPositionFactor.cpp.

21  : IFormFactorDecorator(ff), m_position(position)
22 {
23  setName("FormFactorDecoratorPositionFactor");
24 }
IFormFactorDecorator(const IFormFactor &ff)
void setName(const std::string &name)

References IParametricComponent::setName().

Referenced by clone().

Here is the call graph for this function:

Member Function Documentation

◆ accept()

void FormFactorDecoratorPositionFactor::accept ( INodeVisitor visitor) const
inlinefinalvirtual

Calls the INodeVisitor's visit method.

Implements INode.

Definition at line 37 of file FormFactorDecoratorPositionFactor.h.

37 { visitor->visit(this); }
virtual void visit(const BasicLattice2D *)
Definition: INodeVisitor.h:151

◆ bottomZ()

double FormFactorDecoratorPositionFactor::bottomZ ( const IRotation rotation) const
finalvirtual

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

Implements IFormFactor.

Definition at line 26 of file FormFactorDecoratorPositionFactor.cpp.

27 {
28  kvector_t rotated_translation = rotation.transformed(m_position);
29  return m_ff->bottomZ(rotation) + rotated_translation.z();
30 }
T z() const
Returns z-component in cartesian coordinate system.
Definition: BasicVector3D.h:67
virtual double bottomZ(const IRotation &rotation) const =0
Returns the z-coordinate of the lowest point in this shape after a given rotation.
kvector_t transformed(const kvector_t &v) const
Definition: Rotations.cpp:53

References IFormFactor::bottomZ(), IFormFactorDecorator::m_ff, m_position, IRotation::transformed(), and BasicVector3D< T >::z().

Here is the call graph for this function:

◆ canSliceAnalytically()

bool IFormFactor::canSliceAnalytically ( const IRotation rot) const
protectedvirtualinherited

Checks if slicing has a fast analytical solution.

Reimplemented in IBornFF, FormFactorFullSphere, and FormFactorDot.

Definition at line 89 of file IFormFactor.cpp.

90 {
91  return false;
92 }

Referenced by IFormFactor::createSlicedFormFactor().

◆ clone()

FormFactorDecoratorPositionFactor* FormFactorDecoratorPositionFactor::clone ( ) const
inlinefinalvirtual

Returns a clone of this ISampleNode object.

Implements IFormFactorDecorator.

Definition at line 32 of file FormFactorDecoratorPositionFactor.h.

33  {
35  }
FormFactorDecoratorPositionFactor(const IFormFactor &ff, const kvector_t &position)

References FormFactorDecoratorPositionFactor(), IFormFactorDecorator::m_ff, and m_position.

Here is the call graph for this function:

◆ containedMaterials()

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

Returns set of unique materials contained in this ISampleNode.

Definition at line 26 of file ISampleNode.cpp.

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

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

Referenced by MultiLayerUtils::ContainsCompatibleMaterials(), SampleToPython::initLabels(), and ISampleNode::isMagnetic().

Here is the call graph for this function:

◆ 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 94 of file INode.cpp.

95 {
96  if (node->parent() != this)
97  return -1;
98 
99  int result(-1), count(0);
100  for (auto child : getChildren()) {
101 
102  if (child == nullptr)
103  throw std::runtime_error("INode::copyNumber() -> Error. Nullptr as child.");
104 
105  if (child == node)
106  result = count;
107 
108  if (child->getName() == node->getName())
109  ++count;
110  }
111 
112  return count > 1 ? result : -1;
113 }
const INode * parent() const
Definition: INode.cpp:84
const std::string & getName() const

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

Referenced by INode::displayName().

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

Definition at line 126 of file INode.cpp.

127 {
128  std::unique_ptr<ParameterPool> result(new ParameterPool);
129 
130  for (const INode* child : progeny()) {
131  const std::string path = NodeUtils::nodePath(child, parent()) + "/";
132  child->parameterPool()->copyToExternalPool(path, result.get());
133  }
134 
135  return result.release();
136 }
Base class for tree-like structures containing parameterized objects.
Definition: INode.h:49
std::vector< const INode * > progeny() const
Returns a vector of all descendants.
Definition: INode.cpp:68
Container with parameters for IParametricComponent object.
Definition: ParameterPool.h:29
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:89

References NodeUtils::nodePath(), INode::parent(), and INode::progeny().

Referenced by ISimulation::runSimulation(), DepthProbeSimulation::validateParametrization(), OffSpecularSimulation::validateParametrization(), and SpecularSimulation::validateParametrization().

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 58 of file IFormFactor.cpp.

60 {
61  if (shapeIsContainedInLimits(*this, limits, rot, translation))
62  return createTransformedFormFactor(*this, rot, translation);
63  if (shapeOutsideLimits(*this, limits, rot, translation))
64  return nullptr;
65  if (canSliceAnalytically(rot))
66  return sliceFormFactor(limits, rot, translation);
67  throw std::runtime_error(getName()
68  + "::createSlicedFormFactor error: not supported for "
69  "the given rotation!");
70 }
virtual IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const
Actually slices the form factor or throws an exception.
Definition: IFormFactor.cpp:94
virtual bool canSliceAnalytically(const IRotation &rot) const
Checks if slicing has a fast analytical solution.
Definition: IFormFactor.cpp:89
static IFormFactor * createTransformedFormFactor(const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:99

References IFormFactor::canSliceAnalytically(), IFormFactor::createTransformedFormFactor(), IParametricComponent::getName(), and IFormFactor::sliceFormFactor().

Here is the call graph for this function:

◆ createTransformedFormFactor()

IFormFactor * IFormFactor::createTransformedFormFactor ( const IFormFactor formfactor,
const IRotation rot,
kvector_t  translation 
)
staticprotectedinherited

Definition at line 99 of file IFormFactor.cpp.

101 {
102  std::unique_ptr<IFormFactor> P_fftemp, P_result;
103  if (!rot.isIdentity())
104  P_fftemp = std::make_unique<FormFactorDecoratorRotation>(formfactor, rot);
105  else
106  P_fftemp.reset(formfactor.clone());
107  if (translation != kvector_t())
108  P_result = std::make_unique<FormFactorDecoratorPositionFactor>(*P_fftemp, translation);
109  else
110  std::swap(P_fftemp, P_result);
111  return P_result.release();
112 }
void swap(OutputDataIterator< TValue, TContainer > &left, OutputDataIterator< TValue, TContainer > &right)
make Swappable
BasicVector3D< double > kvector_t
Definition: Vectors3D.h:21
IFormFactor * clone() const override=0
Returns a clone of this ISampleNode object.
virtual bool isIdentity() const
Returns true if rotation matrix is identity matrix (no rotations)
Definition: Rotations.cpp:58

References IFormFactor::clone(), IRotation::isIdentity(), and swap().

Referenced by IFormFactor::createSlicedFormFactor(), 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:

◆ displayName()

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

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

Definition at line 115 of file INode.cpp.

116 {
117  std::string result = getName();
118  if (m_parent) {
119  int index = m_parent->copyNumber(this);
120  if (index >= 0)
121  result = result + std::to_string(index);
122  }
123  return result;
124 }
const INode * m_parent
Definition: INode.h:83
int copyNumber(const INode *node) const
Returns copyNumber of child, which takes into account existence of children with same name.
Definition: INode.cpp:94

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

Referenced by NodeUtils::nodePath().

Here is the call graph for this function:

◆ evaluate()

complex_t FormFactorDecoratorPositionFactor::evaluate ( const WavevectorInfo wavevectors) const
finalvirtual

Returns scattering amplitude for complex wavevectors ki, kf.

Implements IFormFactor.

Definition at line 38 of file FormFactorDecoratorPositionFactor.cpp.

39 {
40  return getPositionFactor(wavevectors) * m_ff->evaluate(wavevectors);
41 }
complex_t getPositionFactor(const WavevectorInfo &wavevectors) const
virtual complex_t evaluate(const WavevectorInfo &wavevectors) const =0
Returns scattering amplitude for complex wavevectors ki, kf.

References IFormFactor::evaluate(), getPositionFactor(), and IFormFactorDecorator::m_ff.

Here is the call graph for this function:

◆ evaluatePol()

Eigen::Matrix2cd FormFactorDecoratorPositionFactor::evaluatePol ( const WavevectorInfo wavevectors) const
finalvirtual

Returns scattering amplitude for matrix interactions.

Reimplemented from IFormFactor.

Definition at line 44 of file FormFactorDecoratorPositionFactor.cpp.

45 {
46  return getPositionFactor(wavevectors) * m_ff->evaluatePol(wavevectors);
47 }
virtual Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const
Returns scattering amplitude for matrix interactions.
Definition: IFormFactor.cpp:72

References IFormFactor::evaluatePol(), getPositionFactor(), and IFormFactorDecorator::m_ff.

Here is the call graph for this function:

◆ getChildren()

◆ getFormFactor()

const IFormFactor* IFormFactorDecorator::getFormFactor ( ) const
inlineinherited

Definition at line 48 of file IFormFactorDecorator.h.

48 { return m_ff; }

References IFormFactorDecorator::m_ff.

◆ getName()

◆ getPositionFactor()

complex_t FormFactorDecoratorPositionFactor::getPositionFactor ( const WavevectorInfo wavevectors) const
private

Definition at line 50 of file FormFactorDecoratorPositionFactor.cpp.

51 {
52  cvector_t q = wavevectors.getQ();
53  return exp_I(m_position.dot(q));
54 }
complex_t exp_I(complex_t z)
Returns exp(I*z), where I is the imaginary unit.
Definition: Complex.h:30
auto dot(const BasicVector3D< U > &v) const
Returns dot product of vectors (antilinear in the first [=self] argument).
cvector_t getQ() const

References BasicVector3D< T >::dot(), exp_I(), WavevectorInfo::getQ(), and m_position.

Referenced by evaluate(), and evaluatePol().

Here is the call graph for this function:

◆ isMagnetic()

bool ISampleNode::isMagnetic ( ) const
inherited

Returns true if there is any magnetic material in this ISampleNode.

Definition at line 40 of file ISampleNode.cpp.

41 {
42  const auto materials = containedMaterials();
43  return std::any_of(materials.cbegin(), materials.cend(),
44  [](const Material* mat) { return mat->isMagneticMaterial(); });
45 }
std::vector< const Material * > containedMaterials() const
Returns set of unique materials contained in this ISampleNode.
Definition: ISampleNode.cpp:26

References ISampleNode::containedMaterials().

Referenced by ProcessedSample::initLayouts().

Here is the call graph for this function:

◆ material()

virtual const Material* ISampleNode::material ( ) const
inlinevirtualinherited

◆ onChange()

◆ parameter()

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

Returns parameter with given 'name'.

Definition at line 86 of file IParametricComponent.cpp.

87 {
88  return m_pool->parameter(name);
89 }
std::unique_ptr< ParameterPool > m_pool
parameter pool (kind of pointer-to-implementation)
QString const & name(EShape k)
Definition: particles.cpp:21

References IParametricComponent::m_pool, and RealSpace::Particles::name().

Referenced by DepthProbeSimulation::initialize(), SpecularSimulation::initialize(), Lattice3D::initialize(), IParticle::registerAbundance(), ParticleLayout::registerParticleDensity(), IParticle::registerPosition(), Layer::registerThickness(), Lattice2D::setRotationEnabled(), and DistributionLogNormal::setUnits().

Here is the call graph for this function:

◆ parameterPool()

ParameterPool* IParametricComponent::parameterPool ( ) const
inlineinherited

Returns pointer to the parameter pool.

Definition at line 39 of file IParametricComponent.h.

39 { return m_pool.get(); } // has non-const usages!

References IParametricComponent::m_pool.

Referenced by INode::INode(), IParametricComponent::IParametricComponent(), pyfmt2::argumentList(), SampleBuilderNode::borrow_builder_parameters(), SampleBuilderNode::reset(), and IDistribution1D::setUnits().

◆ parametersToString()

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

Returns multiline string representing available parameters.

Definition at line 43 of file IParametricComponent.cpp.

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

References IParametricComponent::createParameterTree().

Here is the call graph for this function:

◆ parent() [1/2]

INode * INode::parent ( )
inherited

Definition at line 89 of file INode.cpp.

90 {
91  return const_cast<INode*>(m_parent);
92 }

References INode::m_parent.

◆ parent() [2/2]

◆ progeny()

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

Returns a vector of all descendants.

Definition at line 68 of file INode.cpp.

69 {
70  std::vector<const INode*> result;
71  result.push_back(this);
72  for (const auto* child : getChildren()) {
73  for (const auto* p : child->progeny())
74  result.push_back(p);
75  }
76  return result;
77 }

References INode::getChildren().

Referenced by INode::createParameterTree(), and ParticleDistribution::generateParticles().

Here is the call graph for this function:

◆ radialExtension()

double IFormFactorDecorator::radialExtension ( ) const
inlineoverridevirtualinherited

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 42 of file IFormFactorDecorator.h.

42 { return m_ff->radialExtension(); }
virtual double radialExtension() const =0
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.

References IFormFactorDecorator::m_ff, and IFormFactor::radialExtension().

Here is the call graph for this function:

◆ registerChild()

void INode::registerChild ( INode node)
inherited

Definition at line 57 of file INode.cpp.

58 {
59  ASSERT(node);
60  node->setParent(this);
61 }
#define ASSERT(condition)
Definition: Assert.h:31
virtual void setParent(const INode *newParent)
Definition: INode.cpp:79

References ASSERT, and INode::setParent().

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

Here is the call graph for this function:

◆ registerParameter()

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

Definition at line 51 of file IParametricComponent.cpp.

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

References IParametricComponent::getName(), IParametricComponent::m_pool, RealSpace::Particles::name(), and IParametricComponent::onChange().

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

Here is the call graph for this function:

◆ registerVector()

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

Definition at line 57 of file IParametricComponent.cpp.

59 {
60  registerParameter(XComponentName(base_name), &((*p_vec)[0])).setUnit(units);
61  registerParameter(YComponentName(base_name), &((*p_vec)[1])).setUnit(units);
62  registerParameter(ZComponentName(base_name), &((*p_vec)[2])).setUnit(units);
63 }
static std::string XComponentName(const std::string &base_name)
static std::string ZComponentName(const std::string &base_name)
RealParameter & registerParameter(const std::string &name, double *parpointer)
static std::string YComponentName(const std::string &base_name)
RealParameter & setUnit(const std::string &name)
MVVM_MODEL_EXPORT std::string base_name(const std::string &path)
Provide the filename of a file path.
Definition: fileutils.cpp:78

References ModelView::Utils::base_name(), IParametricComponent::registerParameter(), RealParameter::setUnit(), IParametricComponent::XComponentName(), IParametricComponent::YComponentName(), and IParametricComponent::ZComponentName().

Referenced by Beam::Beam(), DetectionProperties::DetectionProperties(), InterferenceFunctionTwin::InterferenceFunctionTwin(), MultiLayer::MultiLayer(), Lattice3D::initialize(), and IParticle::registerPosition().

Here is the call graph for this function:

◆ removeParameter()

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

Definition at line 91 of file IParametricComponent.cpp.

92 {
93  m_pool->removeParameter(name);
94 }

References IParametricComponent::m_pool, and RealSpace::Particles::name().

Referenced by IParticle::registerAbundance(), ParticleLayout::registerParticleDensity(), Layer::registerThickness(), IParametricComponent::removeVector(), and Lattice2D::setRotationEnabled().

Here is the call graph for this function:

◆ removeVector()

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

Definition at line 96 of file IParametricComponent.cpp.

References ModelView::Utils::base_name(), IParametricComponent::removeParameter(), IParametricComponent::XComponentName(), IParametricComponent::YComponentName(), and IParametricComponent::ZComponentName().

Referenced by IParticle::registerPosition().

Here is the call graph for this function:

◆ setAmbientMaterial()

void IFormFactorDecorator::setAmbientMaterial ( const Material )
inlineoverridevirtualinherited

Passes the material in which this particle is embedded.

Implements IFormFactor.

Definition at line 35 of file IFormFactorDecorator.h.

36  {
38  }
virtual void setAmbientMaterial(const Material &)=0
Passes the material in which this particle is embedded.

References IFormFactorDecorator::m_ff, ISampleNode::material(), and IFormFactor::setAmbientMaterial().

Here is the call graph for this function:

◆ setName()

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

Definition at line 69 of file IParametricComponent.h.

69 { m_name = name; }

References IParametricComponent::m_name, and RealSpace::Particles::name().

Referenced by BasicLattice2D::BasicLattice2D(), Beam::Beam(), ConvolutionDetectorResolution::ConvolutionDetectorResolution(), Crystal::Crystal(), DetectionProperties::DetectionProperties(), DistributionHandler::DistributionHandler(), FormFactorCoreShell::FormFactorCoreShell(), FormFactorCrystal::FormFactorCrystal(), FormFactorDecoratorMaterial::FormFactorDecoratorMaterial(), FormFactorDecoratorPositionFactor(), FormFactorDecoratorRotation::FormFactorDecoratorRotation(), FormFactorWeighted::FormFactorWeighted(), HexagonalLattice2D::HexagonalLattice2D(), IDetector::IDetector(), 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(), Lattice3D::Lattice3D(), Layer::Layer(), LayerInterface::LayerInterface(), LayerRoughness::LayerRoughness(), MultiLayer::MultiLayer(), ParticleCoreShell::ParticleCoreShell(), ParticleDistribution::ParticleDistribution(), ParticleLayout::ParticleLayout(), RectangularDetector::RectangularDetector(), ResolutionFunction2DGaussian::ResolutionFunction2DGaussian(), SampleBuilderNode::SampleBuilderNode(), SphericalDetector::SphericalDetector(), SquareLattice2D::SquareLattice2D(), Layer::clone(), LayersWithAbsorptionBuilder::createSampleByIndex(), Basic2DParaCrystalBuilder::createSampleByIndex(), ParticleInVacuumBuilder::createSampleByIndex(), SimpleMagneticRotationBuilder::createSampleByIndex(), DepthProbeSimulation::initialize(), GISASSimulation::initialize(), OffSpecularSimulation::initialize(), SpecularSimulation::initialize(), SpecularDetector1D::initialize(), MesoCrystal::initialize(), Particle::initialize(), ParticleComposition::initialize(), Beam::operator=(), SampleBuilderNode::operator=(), SampleBuilderNode::reset(), and SampleBuilderNode::setSBN().

Here is the call graph for this function:

◆ setParameterValue()

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

Definition at line 65 of file IParametricComponent.cpp.

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

References IParametricComponent::createParameterTree(), IParametricComponent::m_pool, RealSpace::Particles::name(), and ParameterPool::setMatchedParametersValue().

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

Here is the call graph for this function:

◆ setParent()

void INode::setParent ( const INode newParent)
virtualinherited

Reimplemented in SampleProvider.

Definition at line 79 of file INode.cpp.

80 {
81  m_parent = newParent;
82 }

References INode::m_parent.

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

◆ setSpecularInfo()

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

Sets reflection/transmission info.

Definition at line 84 of file IFormFactor.cpp.

86 {
87 }

◆ setVectorValue()

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

Definition at line 78 of file IParametricComponent.cpp.

79 {
83 }
T y() const
Returns y-component in cartesian coordinate system.
Definition: BasicVector3D.h:65
T x() const
Returns x-component in cartesian coordinate system.
Definition: BasicVector3D.h:63
void setParameterValue(const std::string &name, double value)

References ModelView::Utils::base_name(), IParametricComponent::setParameterValue(), BasicVector3D< T >::x(), IParametricComponent::XComponentName(), BasicVector3D< T >::y(), IParametricComponent::YComponentName(), BasicVector3D< T >::z(), and IParametricComponent::ZComponentName().

Here is the call graph for this function:

◆ 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 94 of file IFormFactor.cpp.

95 {
96  throw std::runtime_error(getName() + "::sliceFormFactor error: not implemented!");
97 }

References IParametricComponent::getName().

Referenced by IFormFactor::createSlicedFormFactor().

Here is the call graph for this function:

◆ topZ()

double FormFactorDecoratorPositionFactor::topZ ( const IRotation rotation) const
finalvirtual

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

Implements IFormFactor.

Definition at line 32 of file FormFactorDecoratorPositionFactor.cpp.

33 {
34  kvector_t rotated_translation = rotation.transformed(m_position);
35  return m_ff->topZ(rotation) + rotated_translation.z();
36 }
virtual double topZ(const IRotation &rotation) const =0
Returns the z-coordinate of the lowest point in this shape after a given rotation.

References IFormFactorDecorator::m_ff, m_position, IFormFactor::topZ(), IRotation::transformed(), and BasicVector3D< T >::z().

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 52 of file INode.cpp.

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

References NodeUtils::nodeToString().

Here is the call graph for this function:

◆ volume()

double IFormFactorDecorator::volume ( ) const
inlineoverridevirtualinherited

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

Reimplemented from IFormFactor.

Definition at line 40 of file IFormFactorDecorator.h.

40 { return m_ff->volume(); }
virtual double volume() const
Returns the total volume of the particle of this form factor's shape.
Definition: IFormFactor.cpp:78

References IFormFactorDecorator::m_ff, and IFormFactor::volume().

Here is the call graph for this function:

◆ XComponentName()

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

Definition at line 103 of file IParametricComponent.cpp.

104 {
105  return base_name + "X";
106 }

References ModelView::Utils::base_name().

Referenced by Lattice3D::initialize(), IParticle::registerPosition(), IParametricComponent::registerVector(), IParametricComponent::removeVector(), IParametricComponent::setVectorValue(), and VectorParameterTranslator::translate().

Here is the call graph for this function:

◆ YComponentName()

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

Definition at line 108 of file IParametricComponent.cpp.

109 {
110  return base_name + "Y";
111 }

References ModelView::Utils::base_name().

Referenced by IParametricComponent::registerVector(), IParametricComponent::removeVector(), IParametricComponent::setVectorValue(), and VectorParameterTranslator::translate().

Here is the call graph for this function:

◆ ZComponentName()

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

Definition at line 113 of file IParametricComponent.cpp.

114 {
115  return base_name + "Z";
116 }

References ModelView::Utils::base_name().

Referenced by IParametricComponent::registerVector(), IParametricComponent::removeVector(), IParametricComponent::setVectorValue(), and VectorParameterTranslator::translate().

Here is the call graph for this function:

Member Data Documentation

◆ m_ff

◆ m_name

std::string IParametricComponent::m_name
privateinherited

◆ m_NP

const size_t INode::m_NP
protectedinherited

Definition at line 88 of file INode.h.

Referenced by INode::INode().

◆ m_P

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

Definition at line 89 of file INode.h.

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

◆ m_parent

const INode* INode::m_parent {nullptr}
privateinherited

Definition at line 83 of file INode.h.

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

◆ m_pool

◆ m_position

kvector_t FormFactorDecoratorPositionFactor::m_position
private

Definition at line 51 of file FormFactorDecoratorPositionFactor.h.

Referenced by bottomZ(), clone(), getPositionFactor(), and topZ().


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