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

Public Member Functions

 InterferenceFunction2DParaCrystal (const Lattice2D &lattice, double damping_length, double domain_size_1, double domain_size_2)
 
 InterferenceFunction2DParaCrystal (double length_1, double length_2, double alpha, double xi, double damping_length)
 
 ~InterferenceFunction2DParaCrystal () final
 
InterferenceFunction2DParaCrystalclone () const override final
 
void accept (INodeVisitor *visitor) const override final
 
void setDomainSizes (double size_1, double size_2)
 
void setProbabilityDistributions (const IFTDistribution2D &pdf_1, const IFTDistribution2D &pdf_2)
 
void setDampingLength (double damping_length)
 
std::vector< double > domainSizes () const
 
void setIntegrationOverXi (bool integrate_xi)
 
bool integrationOverXi () const
 
double dampingLength () const
 
const Lattice2Dlattice () const
 
double getParticleDensity () const override final
 
std::vector< const INode * > getChildren () const override final
 
const IFTDistribution2Dpdf1 () const
 
const IFTDistribution2Dpdf2 () const
 
virtual double evaluate (const kvector_t q, double outer_iff=1.0) const
 
void setPositionVariance (double var)
 
double positionVariance () const
 
virtual bool supportsMultilayer () const
 
double DWfactor (kvector_t q) const
 
virtual const Materialmaterial () const
 
std::vector< const Material * > containedMaterials () const
 
virtual void transferToCPP ()
 
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 InterferenceFunction2DParaCrystalcreateSquare (double lattice_length, double damping_length, double domain_size_1, double domain_size_2)
 
static InterferenceFunction2DParaCrystalcreateHexagonal (double lattice_length, double damping_length, double domain_size_1, double domain_size_2)
 
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

double iff_no_inner (const kvector_t q, double outer_iff) const
 

Protected Attributes

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

Private Member Functions

double iff_without_dw (const kvector_t q) const override final
 
void setLattice (const Lattice2D &lattice)
 
double interferenceForXi (double xi) const
 
double interference1D (double qx, double qy, double xi, size_t index) const
 
complex_t FTPDF (double qx, double qy, double xi, size_t index) const
 
void transformToPrincipalAxes (double qx, double qy, double gamma, double delta, double &q_pa_1, double &q_pa_2) const
 

Private Attributes

bool m_integrate_xi
 
std::unique_ptr< IFTDistribution2Dm_pdf1
 
std::unique_ptr< IFTDistribution2Dm_pdf2
 
std::unique_ptr< Lattice2Dm_lattice
 
double m_damping_length
 
double m_domain_sizes [2]
 
double m_qx
 
double m_qy
 
const INodem_parent {nullptr}
 
std::string m_name
 
std::unique_ptr< ParameterPoolm_pool
 

Detailed Description

Interference function of a 2D paracrystal.

Definition at line 29 of file InterferenceFunction2DParaCrystal.h.

Constructor & Destructor Documentation

◆ InterferenceFunction2DParaCrystal() [1/2]

InterferenceFunction2DParaCrystal::InterferenceFunction2DParaCrystal ( const Lattice2D lattice,
double  damping_length,
double  domain_size_1,
double  domain_size_2 
)

Definition at line 22 of file InterferenceFunction2DParaCrystal.cpp.

26  : IInterferenceFunction(0), m_integrate_xi(false), m_damping_length(damping_length)
27 {
28  setName("Interference2DParaCrystal");
30  setDomainSizes(domain_size_1, domain_size_2);
31  registerParameter("DampingLength", &m_damping_length).setUnit("nm").setNonnegative();
32  registerParameter("DomainSize1", &m_domain_sizes[0]).setUnit("nm").setNonnegative();
33  registerParameter("DomainSize2", &m_domain_sizes[1]).setUnit("nm").setNonnegative();
34 }
IInterferenceFunction(const NodeMeta &meta, const std::vector< double > &PValues)
RealParameter & registerParameter(const std::string &name, double *parpointer)
void setName(const std::string &name)
double m_domain_sizes[2]
Coherence domain sizes.
bool m_integrate_xi
Integrate over the orientation xi.
double m_damping_length
Damping length for removing delta function singularity at q=0.
void setDomainSizes(double size_1, double size_2)
Sets the sizes of coherence domains.
RealParameter & setNonnegative()
RealParameter & setUnit(const std::string &name)

References lattice(), m_damping_length, m_domain_sizes, IParameterized::registerParameter(), setDomainSizes(), setLattice(), IParameterized::setName(), RealParameter::setNonnegative(), and RealParameter::setUnit().

Referenced by clone(), createHexagonal(), and createSquare().

Here is the call graph for this function:

◆ InterferenceFunction2DParaCrystal() [2/2]

InterferenceFunction2DParaCrystal::InterferenceFunction2DParaCrystal ( double  length_1,
double  length_2,
double  alpha,
double  xi,
double  damping_length 
)

Constructor of interference function of two-dimensional paracrystal.

Parameters
length_1length of first lattice vector in nanometers
length_2length of second lattice vector in nanometers
alphaangle between lattice vectors in radians
xirotation of lattice with respect to x-axis (beam direction) in radians
damping_lengththe damping (coherence) length of the paracrystal in nanometers

Definition at line 43 of file InterferenceFunction2DParaCrystal.cpp.

47  : InterferenceFunction2DParaCrystal(BasicLattice(length_1, length_2, alpha, xi), damping_length,
48  0, 0)
49 {
50 }
InterferenceFunction2DParaCrystal(const Lattice2D &lattice, double damping_length, double domain_size_1, double domain_size_2)

◆ ~InterferenceFunction2DParaCrystal()

InterferenceFunction2DParaCrystal::~InterferenceFunction2DParaCrystal ( )
finaldefault

Member Function Documentation

◆ clone()

InterferenceFunction2DParaCrystal * InterferenceFunction2DParaCrystal::clone ( ) const
finaloverridevirtual

Returns a clone of this ISample object.

Implements IInterferenceFunction.

Definition at line 54 of file InterferenceFunction2DParaCrystal.cpp.

55 {
58  ret->setPositionVariance(m_position_var);
59  if (m_pdf1 && m_pdf2)
60  ret->setProbabilityDistributions(*m_pdf1, *m_pdf2);
61  ret->setIntegrationOverXi(m_integrate_xi);
62  return ret;
63 }
std::unique_ptr< IFTDistribution2D > m_pdf1
std::unique_ptr< IFTDistribution2D > m_pdf2

References InterferenceFunction2DParaCrystal(), m_damping_length, m_domain_sizes, m_integrate_xi, m_lattice, m_pdf1, m_pdf2, and IInterferenceFunction::m_position_var.

Here is the call graph for this function:

◆ accept()

void InterferenceFunction2DParaCrystal::accept ( INodeVisitor visitor) const
inlinefinaloverridevirtual

Calls the INodeVisitor's visit method.

Implements INode.

Definition at line 42 of file InterferenceFunction2DParaCrystal.h.

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

◆ createSquare()

InterferenceFunction2DParaCrystal * InterferenceFunction2DParaCrystal::createSquare ( double  lattice_length,
double  damping_length,
double  domain_size_1,
double  domain_size_2 
)
static

Creates square lattice.

Parameters
lattice_lengthlength of first and second lattice vectors in nanometers
damping_lengththe damping (coherence) length of the paracrystal in nanometers
domain_size_1size of the coherent domain along the first basis vector in nanometers
domain_size_2size of the coherent domain along the second basis vector in nanometers

Definition at line 121 of file InterferenceFunction2DParaCrystal.cpp.

123 {
124  auto result = new InterferenceFunction2DParaCrystal(
125  SquareLattice(lattice_length), damping_length, domain_size_1, domain_size_2);
126  result->setIntegrationOverXi(true);
127  return result;
128 }

References InterferenceFunction2DParaCrystal().

Referenced by RectParaCrystalBuilder::buildSample().

Here is the call graph for this function:

◆ createHexagonal()

InterferenceFunction2DParaCrystal * InterferenceFunction2DParaCrystal::createHexagonal ( double  lattice_length,
double  damping_length,
double  domain_size_1,
double  domain_size_2 
)
static

Creates hexagonal lattice.

Parameters
lattice_lengthlength of first and second lattice vectors in nanometers
damping_lengththe damping (coherence) length of the paracrystal in nanometers
domain_size_1size of the coherent domain along the first basis vector in nanometers
domain_size_2size of the coherent domain along the second basis vector in nanometers

Definition at line 137 of file InterferenceFunction2DParaCrystal.cpp.

139 {
140  auto result = new InterferenceFunction2DParaCrystal(
141  HexagonalLattice(lattice_length, 0.), damping_length, domain_size_1, domain_size_2);
142  result->setIntegrationOverXi(true);
143  return result;
144 }

References InterferenceFunction2DParaCrystal().

Referenced by HexParaCrystalBuilder::buildSample().

Here is the call graph for this function:

◆ setDomainSizes()

void InterferenceFunction2DParaCrystal::setDomainSizes ( double  size_1,
double  size_2 
)

Sets the sizes of coherence domains.

Parameters
size_1coherence domain size along the first basis vector in nanometers
size_2coherence domain size along the second basis vector in nanometers

Definition at line 150 of file InterferenceFunction2DParaCrystal.cpp.

151 {
152  m_domain_sizes[0] = size_1;
153  m_domain_sizes[1] = size_2;
154 }

References m_domain_sizes.

Referenced by Basic2DParaCrystalBuilder::buildSample(), RectParaCrystalBuilder::buildSample(), and InterferenceFunction2DParaCrystal().

◆ setProbabilityDistributions()

void InterferenceFunction2DParaCrystal::setProbabilityDistributions ( const IFTDistribution2D pdf_1,
const IFTDistribution2D pdf_2 
)

Sets the probability distributions (Fourier transformed) for the two lattice directions.

Parameters
pdf_1probability distribution in first lattice direction
pdf_2probability distribution in second lattice direction

Definition at line 69 of file InterferenceFunction2DParaCrystal.cpp.

71 {
72  m_pdf1.reset(pdf_1.clone());
73  registerChild(m_pdf1.get());
74  m_pdf2.reset(pdf_2.clone());
75  registerChild(m_pdf2.get());
76 }
IFTDistribution2D * clone() const =0
void registerChild(INode *node)
Definition: INode.cpp:58

References IFTDistribution2D::clone(), m_pdf1, m_pdf2, and INode::registerChild().

Referenced by Basic2DParaCrystalBuilder::buildSample().

Here is the call graph for this function:

◆ setDampingLength()

void InterferenceFunction2DParaCrystal::setDampingLength ( double  damping_length)

Sets the damping length.

Parameters
damping_lengththe damping (coherence) length of the paracrystal in nanometers

Definition at line 81 of file InterferenceFunction2DParaCrystal.cpp.

82 {
83  m_damping_length = damping_length;
84 }

References m_damping_length.

◆ domainSizes()

std::vector< double > InterferenceFunction2DParaCrystal::domainSizes ( ) const

Definition at line 233 of file InterferenceFunction2DParaCrystal.cpp.

234 {
235  return {m_domain_sizes[0], m_domain_sizes[1]};
236 }

References m_domain_sizes.

◆ setIntegrationOverXi()

void InterferenceFunction2DParaCrystal::setIntegrationOverXi ( bool  integrate_xi)

Enables/disables averaging over the lattice rotation angle.

Parameters
integrate_xiintegration flag

Definition at line 241 of file InterferenceFunction2DParaCrystal.cpp.

242 {
243  m_integrate_xi = integrate_xi;
244  m_lattice->setRotationEnabled(!m_integrate_xi); // deregister Xi in the case of integration
245 }

References m_integrate_xi, and m_lattice.

◆ integrationOverXi()

bool InterferenceFunction2DParaCrystal::integrationOverXi ( ) const
inline

Definition at line 63 of file InterferenceFunction2DParaCrystal.h.

63 { return m_integrate_xi; }

References m_integrate_xi.

◆ dampingLength()

double InterferenceFunction2DParaCrystal::dampingLength ( ) const
inline

Definition at line 64 of file InterferenceFunction2DParaCrystal.h.

64 { return m_damping_length; }

References m_damping_length.

◆ lattice()

const Lattice2D & InterferenceFunction2DParaCrystal::lattice ( ) const

Definition at line 247 of file InterferenceFunction2DParaCrystal.cpp.

248 {
249  if (!m_lattice)
250  throw std::runtime_error("InterferenceFunction2DParaCrystal::lattice() -> Error. "
251  "No lattice defined.");
252  return *m_lattice;
253 }

References m_lattice.

Referenced by InterferenceFunction2DParaCrystal(), and setLattice().

◆ getParticleDensity()

double InterferenceFunction2DParaCrystal::getParticleDensity ( ) const
finaloverridevirtual

If defined by this interference function's parameters, returns the particle density (per area).

Otherwise, returns zero or a user-defined value

Reimplemented from IInterferenceFunction.

Definition at line 86 of file InterferenceFunction2DParaCrystal.cpp.

87 {
88  double area = m_lattice->unitCellArea();
89  return area == 0.0 ? 0.0 : 1.0 / area;
90 }

References m_lattice.

◆ getChildren()

std::vector< const INode * > InterferenceFunction2DParaCrystal::getChildren ( ) const
finaloverridevirtual

Returns a vector of children (const).

Reimplemented from INode.

Definition at line 92 of file InterferenceFunction2DParaCrystal.cpp.

93 {
94  return std::vector<const INode*>() << m_pdf1 << m_pdf2 << m_lattice;
95 }

References m_lattice, m_pdf1, and m_pdf2.

◆ pdf1()

const IFTDistribution2D* InterferenceFunction2DParaCrystal::pdf1 ( ) const
inline

Definition at line 72 of file InterferenceFunction2DParaCrystal.h.

72 { return m_pdf1.get(); }

References m_pdf1.

◆ pdf2()

const IFTDistribution2D* InterferenceFunction2DParaCrystal::pdf2 ( ) const
inline

Definition at line 74 of file InterferenceFunction2DParaCrystal.h.

74 { return m_pdf2.get(); }

References m_pdf2.

◆ iff_without_dw()

double InterferenceFunction2DParaCrystal::iff_without_dw ( const kvector_t  q) const
finaloverrideprivatevirtual

Calculates the structure factor without Debye-Waller factor.

Implements IInterferenceFunction.

Definition at line 97 of file InterferenceFunction2DParaCrystal.cpp.

98 {
99  m_qx = q.x();
100  m_qy = q.y();
101  if (!m_integrate_xi)
102  return interferenceForXi(m_lattice->rotationAngle());
103  return RealIntegrator().integrate([&](double xi) -> double { return interferenceForXi(xi); },
104  0.0, M_TWOPI)
105  / M_TWOPI;
106 }
#define M_TWOPI
Definition: MathConstants.h:49
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
double interferenceForXi(double xi) const
Returns interference function for fixed angle xi.
To integrate a real function of a real variable.
Definition: Integrator.h:24
double integrate(const std::function< double(double)> &f, double lmin, double lmax)
Definition: Integrator.cpp:27

References RealIntegrator::integrate(), interferenceForXi(), m_integrate_xi, m_lattice, m_qx, m_qy, M_TWOPI, BasicVector3D< T >::x(), and BasicVector3D< T >::y().

Here is the call graph for this function:

◆ setLattice()

void InterferenceFunction2DParaCrystal::setLattice ( const Lattice2D lattice)
private

Definition at line 108 of file InterferenceFunction2DParaCrystal.cpp.

109 {
110  m_lattice.reset(lattice.clone());
111  registerChild(m_lattice.get());
112 }
virtual Lattice2D * clone() const =0

References Lattice2D::clone(), lattice(), m_lattice, and INode::registerChild().

Referenced by InterferenceFunction2DParaCrystal().

Here is the call graph for this function:

◆ interferenceForXi()

double InterferenceFunction2DParaCrystal::interferenceForXi ( double  xi) const
private

Returns interference function for fixed angle xi.

Definition at line 165 of file InterferenceFunction2DParaCrystal.cpp.

166 {
167  // don't touch order of computation; problems under Windows
168  double rx = interference1D(m_qx, m_qy, xi, 0);
169  double ry = interference1D(m_qx, m_qy, xi + m_lattice->latticeAngle(), 1);
170  return rx * ry;
171 }
double interference1D(double qx, double qy, double xi, size_t index) const
Returns interference function for fixed xi in the dimension determined by the given index.

References interference1D(), m_lattice, m_qx, and m_qy.

Referenced by iff_without_dw().

Here is the call graph for this function:

◆ interference1D()

double InterferenceFunction2DParaCrystal::interference1D ( double  qx,
double  qy,
double  xi,
size_t  index 
) const
private

Returns interference function for fixed xi in the dimension determined by the given index.

Definition at line 174 of file InterferenceFunction2DParaCrystal.cpp.

176 {
177  if (index > 1)
179  "InterferenceFunction2DParaCrystal::"
180  "interference1D() -> Error! Index of interference function "
181  "probability must be < 2");
182  if (!m_pdf1 || !m_pdf2)
184  "InterferenceFunction2DParaCrystal::"
185  "interference1D() -> Error! Probability distributions for "
186  "interference function not properly initialized");
187 
188  double length = index ? m_lattice->length2() : m_lattice->length1();
189  int n = static_cast<int>(std::abs(m_domain_sizes[index] / length));
190  double nd = static_cast<double>(n);
191  complex_t fp = FTPDF(qx, qy, xi, index);
192  if (n < 1)
193  return ((1.0 + fp) / (1.0 - fp)).real();
194  if (std::norm(1.0 - fp) < std::numeric_limits<double>::epsilon())
195  return nd;
196  // for (1-fp)*nd small, take the series expansion to second order in nd*(1-fp)
197  if (std::abs(1.0 - fp) * nd < 2e-4) {
198  complex_t intermediate =
199  (nd - 1.0) / 2.0 + (nd * nd - 1.0) * (fp - 1.0) / 6.0
200  + (nd * nd * nd - 2.0 * nd * nd - nd + 2.0) * (fp - 1.0) * (fp - 1.0) / 24.0;
201  return 1.0 + 2.0 * intermediate.real();
202  }
203  complex_t tmp;
204  if (std::abs(fp) == 0.0
205  || std::log(std::abs(fp)) * nd < std::log(std::numeric_limits<double>::min()))
206  tmp = 0.0;
207  else
208  tmp = std::pow(fp, n);
209  complex_t intermediate = fp / (1.0 - fp) - fp * (1.0 - tmp) / nd / (1.0 - fp) / (1.0 - fp);
210  return 1.0 + 2.0 * intermediate.real();
211 }
std::complex< double > complex_t
Definition: Complex.h:20
complex_t FTPDF(double qx, double qy, double xi, size_t index) const

References FTPDF(), anonymous_namespace{BoxCompositionBuilder.cpp}::length, m_domain_sizes, m_lattice, m_pdf1, and m_pdf2.

Referenced by interferenceForXi().

Here is the call graph for this function:

◆ FTPDF()

complex_t InterferenceFunction2DParaCrystal::FTPDF ( double  qx,
double  qy,
double  xi,
size_t  index 
) const
private

Definition at line 213 of file InterferenceFunction2DParaCrystal.cpp.

215 {
216  double length = (index ? m_lattice->length2() : m_lattice->length1());
217 
218  const IFTDistribution2D* pdf = (index ? m_pdf2.get() : m_pdf1.get());
219  double qa = qx * length * std::cos(xi) + qy * length * std::sin(xi);
220  complex_t phase = exp_I(qa);
221  // transform q to principal axes:
222  double qp1, qp2;
223  double gamma = xi + pdf->gamma();
224  double delta = pdf->delta();
225  transformToPrincipalAxes(qx, qy, gamma, delta, qp1, qp2);
226  double amplitude = pdf->evaluate(qp1, qp2);
227  complex_t result = phase * amplitude;
228  if (m_damping_length != 0.0)
229  result *= std::exp(-length / m_damping_length);
230  return result;
231 }
complex_t exp_I(complex_t z)
Returns exp(I*z), where I is the imaginary unit.
Definition: Complex.h:30
Interface for two-dimensional distributions in Fourier space.
double delta() const
Angle in direct space between X- and Y-axis of distribution.
virtual double evaluate(double qx, double qy) const =0
evaluate Fourier transformed distribution for q in X,Y coordinates the original distribution (in real...
double gamma() const
void transformToPrincipalAxes(double qx, double qy, double gamma, double delta, double &q_pa_1, double &q_pa_2) const

References IFTDistribution2D::delta(), IFTDistribution2D::evaluate(), exp_I(), IFTDistribution2D::gamma(), anonymous_namespace{BoxCompositionBuilder.cpp}::length, m_damping_length, m_lattice, m_pdf1, m_pdf2, and transformToPrincipalAxes().

Referenced by interference1D().

Here is the call graph for this function:

◆ transformToPrincipalAxes()

void InterferenceFunction2DParaCrystal::transformToPrincipalAxes ( double  qx,
double  qy,
double  gamma,
double  delta,
double &  q_pa_1,
double &  q_pa_2 
) const
private

Definition at line 156 of file InterferenceFunction2DParaCrystal.cpp.

159 {
160  q_pa_1 = qx * std::cos(gamma) + qy * std::sin(gamma);
161  q_pa_2 = qx * std::cos(gamma + delta) + qy * std::sin(gamma + delta);
162 }

Referenced by FTPDF().

◆ evaluate()

double IInterferenceFunction::evaluate ( const kvector_t  q,
double  outer_iff = 1.0 
) const
virtualinherited

Evaluates the interference function for a given wavevector transfer.

Reimplemented in InterferenceFunction2DSuperLattice.

Definition at line 35 of file IInterferenceFunction.cpp.

36 {
37  return iff_no_inner(q, outer_iff);
38 }
double iff_no_inner(const kvector_t q, double outer_iff) const
Calculates the structure factor in the absence of extra inner structure.

References IInterferenceFunction::iff_no_inner().

Here is the call graph for this function:

◆ setPositionVariance()

void IInterferenceFunction::setPositionVariance ( double  var)
inherited

Sets the variance of the position for the calculation of the DW factor It is defined as the variance in each relevant dimension.

Definition at line 40 of file IInterferenceFunction.cpp.

41 {
42  if (var < 0.0)
43  throw std::runtime_error("IInterferenceFunction::setPositionVariance: "
44  "variance should be positive.");
45  m_position_var = var;
46 }

References IInterferenceFunction::m_position_var.

Referenced by FiniteSquareLatticeBuilder::buildSample(), and SuperLatticeBuilder::buildSample().

◆ positionVariance()

double IInterferenceFunction::positionVariance ( ) const
inlineinherited

Returns the position variance.

Definition at line 39 of file IInterferenceFunction.h.

39 { return m_position_var; }

References IInterferenceFunction::m_position_var.

Referenced by SampleToPython::defineInterferenceFunctions().

◆ supportsMultilayer()

virtual bool IInterferenceFunction::supportsMultilayer ( ) const
inlinevirtualinherited

Indicates if this interference function can be used with a multilayer (DWBA mode)

Reimplemented in InterferenceFunctionFinite3DLattice, and InterferenceFunction3DLattice.

Definition at line 46 of file IInterferenceFunction.h.

46 { return true; }

Referenced by LayoutStrategyBuilder::createStrategy(), and IInterferenceFunction::DWfactor().

◆ DWfactor()

double IInterferenceFunction::DWfactor ( kvector_t  q) const
inherited

Evaluates the Debye-Waller factor for a given wavevector transfer.

Definition at line 48 of file IInterferenceFunction.cpp.

49 {
50  // remove z component for two dimensional interference functions:
51  if (supportsMultilayer())
52  q.setZ(0.0);
53  return std::exp(-q.mag2() * m_position_var);
54 }
double mag2() const
Returns magnitude squared of the vector.
void setZ(const T &a)
Sets z-component in cartesian coordinate system.
Definition: BasicVector3D.h:75
virtual bool supportsMultilayer() const
Indicates if this interference function can be used with a multilayer (DWBA mode)

References IInterferenceFunction::m_position_var, BasicVector3D< T >::mag2(), BasicVector3D< T >::setZ(), and IInterferenceFunction::supportsMultilayer().

Referenced by IInterferenceFunction::iff_no_inner().

Here is the call graph for this function:

◆ iff_no_inner()

double IInterferenceFunction::iff_no_inner ( const kvector_t  q,
double  outer_iff 
) const
protectedinherited

Calculates the structure factor in the absence of extra inner structure.

Definition at line 56 of file IInterferenceFunction.cpp.

57 {
58  return DWfactor(q) * (iff_without_dw(q) * outer_iff - 1.0) + 1.0;
59 }
double DWfactor(kvector_t q) const
Evaluates the Debye-Waller factor for a given wavevector transfer.
virtual double iff_without_dw(const kvector_t q) const =0
Calculates the structure factor without Debye-Waller factor.

References IInterferenceFunction::DWfactor(), and IInterferenceFunction::iff_without_dw().

Referenced by IInterferenceFunction::evaluate(), and InterferenceFunction2DSuperLattice::interferenceForXi().

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 }
#define ASSERT(condition)
Definition: Assert.h:26
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(), setLattice(), InterferenceFunction2DSuperLattice::setLattice(), InterferenceFunctionFinite2DLattice::setLattice(), InterferenceFunctionRadialParaCrystal::setProbabilityDistribution(), 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
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(), 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(), 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 }
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)

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
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(), 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_integrate_xi

bool InterferenceFunction2DParaCrystal::m_integrate_xi
private

Integrate over the orientation xi.

Definition at line 86 of file InterferenceFunction2DParaCrystal.h.

Referenced by clone(), iff_without_dw(), integrationOverXi(), and setIntegrationOverXi().

◆ m_pdf1

std::unique_ptr<IFTDistribution2D> InterferenceFunction2DParaCrystal::m_pdf1
private

◆ m_pdf2

std::unique_ptr<IFTDistribution2D> InterferenceFunction2DParaCrystal::m_pdf2
private

◆ m_lattice

std::unique_ptr<Lattice2D> InterferenceFunction2DParaCrystal::m_lattice
private

◆ m_damping_length

double InterferenceFunction2DParaCrystal::m_damping_length
private

Damping length for removing delta function singularity at q=0.

Definition at line 89 of file InterferenceFunction2DParaCrystal.h.

Referenced by clone(), dampingLength(), FTPDF(), InterferenceFunction2DParaCrystal(), and setDampingLength().

◆ m_domain_sizes

double InterferenceFunction2DParaCrystal::m_domain_sizes[2]
private

Coherence domain sizes.

Definition at line 90 of file InterferenceFunction2DParaCrystal.h.

Referenced by clone(), domainSizes(), interference1D(), InterferenceFunction2DParaCrystal(), and setDomainSizes().

◆ m_qx

double InterferenceFunction2DParaCrystal::m_qx
mutableprivate

Definition at line 91 of file InterferenceFunction2DParaCrystal.h.

Referenced by iff_without_dw(), and interferenceForXi().

◆ m_qy

double InterferenceFunction2DParaCrystal::m_qy
mutableprivate

Definition at line 92 of file InterferenceFunction2DParaCrystal.h.

Referenced by iff_without_dw(), and interferenceForXi().

◆ m_position_var

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