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

Description

Interference function of radial paracrystal.

Definition at line 26 of file InterferenceRadialParaCrystal.h.

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

Public Member Functions

 InterferenceRadialParaCrystal (double peak_distance, double damping_length)
 Constructor of interference function of radial paracrystal. More...
 
void checkNodeArgs () const
 Raises exception if a parameter value is invalid. More...
 
std::string className () const final
 Returns the class name, to be hard-coded in each leaf class that inherits from INode. More...
 
InterferenceRadialParaCrystalclone () const override
 
double dampingLength () const
 
double domainSize () const
 
double DWfactor (R3 q) const
 structureFactors the Debye-Waller factor for a given wavevector transfer More...
 
complex_t FTPDF (double qpar) const
 
double kappa () const
 
std::vector< const INode * > nodeChildren () const override
 Returns all children. More...
 
std::vector< const INode * > nodeOffspring () const
 Returns all descendants. More...
 
std::vector< ParaMetaparDefs () const final
 Returns the parameter definitions, to be hard-coded in each leaf class. More...
 
virtual double particleDensity () const
 If defined by this interference function's parameters, returns the particle density (per area). Otherwise, returns zero or a user-defined value. More...
 
double peakDistance () const
 
double positionVariance () const
 Returns the position variance. More...
 
double randomSample () const
 
void setDomainSize (double size)
 Sets domain size (finite size corrections). More...
 
void setKappa (double kappa)
 Sets size spacing coupling parameter of the Size Spacing Correlation Approximation. More...
 
void setPositionVariance (double var)
 Sets the variance of the position for the calculation of the DW factor It is defined as the variance in each relevant dimension. More...
 
void setProbabilityDistribution (const IProfile1D &pdf)
 Sets one-dimensional probability distribution. More...
 
virtual double structureFactor (R3 q, double outer_iff=1.0) const
 The interference function for a given wavevector transfer. More...
 
virtual bool supportsMultilayer () const
 Indicates if this interference function can be used with a sample (DWBA mode) More...
 
virtual void transferToCPP ()
 Used for Python overriding of clone (see swig/tweaks.py) More...
 

Protected Member Functions

double iff_no_inner (R3 q, double outer_iff) const
 Calculates the structure factor in the absence of extra inner structure. More...
 

Protected Attributes

std::vector< double > m_P
 
double m_position_var
 

Private Member Functions

double iff_without_dw (R3 q) const override
 Calculates the structure factor without Debye-Waller factor. More...
 
void init_parameters ()
 

Private Attributes

double m_damping_length
 damping length of paracrystal Fourier transformed probability distribution of the nearest particle More...
 
double m_domain_size
 Size of coherence domain. More...
 
double m_kappa
 Size-spacing coupling parameter. More...
 
std::unique_ptr< IProfile1Dm_pdf
 
double m_peak_distance
 the distance to the first neighbor peak More...
 
bool m_use_damping_length
 

Constructor & Destructor Documentation

◆ InterferenceRadialParaCrystal()

InterferenceRadialParaCrystal::InterferenceRadialParaCrystal ( double  peak_distance,
double  damping_length 
)

Constructor of interference function of radial paracrystal.

Parameters
peak_distanceaverage distance to the next neighbor in nanometers
damping_lengththe damping (coherence) length of the paracrystal in nanometers

Definition at line 22 of file InterferenceRadialParaCrystal.cpp.

24  : IInterference(0)
25  , m_peak_distance(peak_distance)
26  , m_damping_length(damping_length)
27  , m_use_damping_length(true)
28  , m_kappa(0.0)
29  , m_domain_size(0.0)
30 {
31  if (m_damping_length == 0.0)
32  m_use_damping_length = false;
33 
36  RealLimits::nonnegative().check("SizeSpaceCoupling", m_kappa);
38 }
IInterference(const std::vector< double > &PValues)
double m_kappa
Size-spacing coupling parameter.
double m_domain_size
Size of coherence domain.
double m_damping_length
damping length of paracrystal Fourier transformed probability distribution of the nearest particle
double m_peak_distance
the distance to the first neighbor peak
void check(const std::string &name, double value) const
Throws if value is outside limits. Parameter 'name' is for exception message.
Definition: RealLimits.cpp:170
static RealLimits nonnegative()
Creates an object which can have only positive values with 0. included.
Definition: RealLimits.cpp:124

References RealLimits::check(), m_damping_length, m_domain_size, m_kappa, m_peak_distance, m_use_damping_length, and RealLimits::nonnegative().

Referenced by clone().

Here is the call graph for this function:

Member Function Documentation

◆ checkNodeArgs()

void INode::checkNodeArgs ( ) const
inherited

Raises exception if a parameter value is invalid.

Definition at line 27 of file INode.cpp.

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

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

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

Here is the call graph for this function:

◆ className()

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

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

Implements INode.

Definition at line 30 of file InterferenceRadialParaCrystal.h.

30 { return "InterferenceRadialParaCrystal"; }

◆ clone()

InterferenceRadialParaCrystal * InterferenceRadialParaCrystal::clone ( ) const
overridevirtual

Implements IInterference.

Definition at line 40 of file InterferenceRadialParaCrystal.cpp.

41 {
43  result->setPositionVariance(m_position_var);
44  if (m_pdf)
45  result->setProbabilityDistribution(*m_pdf);
46  result->setKappa(m_kappa);
47  result->setDomainSize(m_domain_size);
48  return result;
49 }
double m_position_var
Definition: IInterference.h:53
std::unique_ptr< IProfile1D > m_pdf
InterferenceRadialParaCrystal(double peak_distance, double damping_length)
Constructor of interference function of radial paracrystal.

References InterferenceRadialParaCrystal(), m_damping_length, m_domain_size, m_kappa, m_pdf, m_peak_distance, and IInterference::m_position_var.

Here is the call graph for this function:

◆ dampingLength()

double InterferenceRadialParaCrystal::dampingLength ( ) const
inline

Definition at line 49 of file InterferenceRadialParaCrystal.h.

49 { return m_damping_length; }

References m_damping_length.

◆ domainSize()

double InterferenceRadialParaCrystal::domainSize ( ) const
inline

Definition at line 41 of file InterferenceRadialParaCrystal.h.

41 { return m_domain_size; }

References m_domain_size.

◆ DWfactor()

double IInterference::DWfactor ( R3  q) const
inherited

structureFactors the Debye-Waller factor for a given wavevector transfer

Definition at line 48 of file IInterference.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 }
virtual bool supportsMultilayer() const
Indicates if this interference function can be used with a sample (DWBA mode)
Definition: IInterference.h:47

References IInterference::m_position_var, and IInterference::supportsMultilayer().

Referenced by IInterference::iff_no_inner(), and Interference2DSuperLattice::interferenceForXi().

Here is the call graph for this function:

◆ FTPDF()

complex_t InterferenceRadialParaCrystal::FTPDF ( double  qpar) const

Definition at line 70 of file InterferenceRadialParaCrystal.cpp.

71 {
72  complex_t phase = exp_I(qpar * m_peak_distance);
73  double amplitude = m_pdf->standardizedFT(qpar);
74  complex_t result = phase * amplitude;
76  result *= std::exp(-m_peak_distance / m_damping_length);
77  return result;
78 }

References m_damping_length, m_pdf, m_peak_distance, and m_use_damping_length.

Referenced by iff_without_dw().

◆ iff_no_inner()

double IInterference::iff_no_inner ( R3  q,
double  outer_iff 
) const
protectedinherited

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

Definition at line 56 of file IInterference.cpp.

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

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

Referenced by IInterference::structureFactor().

Here is the call graph for this function:

◆ iff_without_dw()

double InterferenceRadialParaCrystal::iff_without_dw ( R3  q) const
overrideprivatevirtual

Calculates the structure factor without Debye-Waller factor.

Implements IInterference.

Definition at line 93 of file InterferenceRadialParaCrystal.cpp.

94 {
95  if (!m_pdf)
96  throw std::runtime_error("InterferenceRadialParaCrystal::"
97  "evaluate() -> Error! Probability distribution for "
98  "interference function not properly initialized");
99  double result = 0.0;
100  double qxr = q.x();
101  double qyr = q.y();
102  double qpar = std::sqrt(qxr * qxr + qyr * qyr);
103  int n = static_cast<int>(std::abs(m_domain_size / m_peak_distance));
104  auto nd = static_cast<double>(n);
105  complex_t fp = FTPDF(qpar);
106  if (n < 1) {
107  if (std::abs(1.0 - fp) < 10. * std::numeric_limits<double>::epsilon())
108  result = m_pdf->qSecondDerivative() / m_peak_distance / m_peak_distance;
109  else
110  result = ((1.0 + fp) / (1.0 - fp)).real();
111  } else {
112  if (std::norm(1.0 - fp) < 10. * std::numeric_limits<double>::epsilon())
113  result = nd;
114  // for (1-fp)*nd small, take the series expansion to second order in nd*(1-fp)
115  else if (std::abs(1.0 - fp) * nd < 2e-4) {
116  complex_t intermediate =
117  (nd - 1.0) / 2.0 + (nd * nd - 1.0) * (fp - 1.0) / 6.0
118  + (nd * nd * nd - 2.0 * nd * nd - nd + 2.0) * (fp - 1.0) * (fp - 1.0) / 24.0;
119  result = 1.0 + 2.0 * intermediate.real();
120  } else {
121  complex_t tmp;
122  if (std::abs(fp) == 0.0
123  || std::log(std::abs(fp)) * nd < std::log(std::numeric_limits<double>::min())) {
124  tmp = 0.0;
125  } else {
126  tmp = std::pow(fp, n);
127  }
128  complex_t intermediate =
129  fp / (1.0 - fp) - fp * (1.0 - tmp) / nd / (1.0 - fp) / (1.0 - fp);
130  result = 1.0 + 2.0 * intermediate.real();
131  }
132  }
133  return result;
134 }

References FTPDF(), m_domain_size, m_pdf, and m_peak_distance.

Here is the call graph for this function:

◆ init_parameters()

void InterferenceRadialParaCrystal::init_parameters ( )
private

◆ kappa()

double InterferenceRadialParaCrystal::kappa ( ) const

Definition at line 57 of file InterferenceRadialParaCrystal.cpp.

58 {
59  return m_kappa;
60 }

References m_kappa.

Referenced by setKappa().

◆ nodeChildren()

std::vector< const INode * > InterferenceRadialParaCrystal::nodeChildren ( ) const
overridevirtual

Returns all children.

Reimplemented from INode.

Definition at line 88 of file InterferenceRadialParaCrystal.cpp.

89 {
90  return std::vector<const INode*>() << m_pdf;
91 }

References m_pdf.

◆ nodeOffspring()

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

Returns all descendants.

Definition at line 61 of file INode.cpp.

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

References INode::nodeChildren().

Here is the call graph for this function:

◆ parDefs()

std::vector<ParaMeta> InterferenceRadialParaCrystal::parDefs ( ) const
inlinefinalvirtual

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

Reimplemented from INode.

Definition at line 31 of file InterferenceRadialParaCrystal.h.

32  {
33  return {{"PeakDistance", "nm", "peak distance", 0, +INF, 0},
34  {"DampingLength", "nm", "damping length", 0, +INF, 0}};
35  }

References INF.

◆ particleDensity()

virtual double IInterference::particleDensity ( ) const
inlinevirtualinherited

If defined by this interference function's parameters, returns the particle density (per area). Otherwise, returns zero or a user-defined value.

Reimplemented in InterferenceHardDisk, InterferenceFinite2DLattice, Interference2DParaCrystal, and Interference2DLattice.

Definition at line 44 of file IInterference.h.

44 { return 0.0; }

◆ peakDistance()

double InterferenceRadialParaCrystal::peakDistance ( ) const
inline

Definition at line 47 of file InterferenceRadialParaCrystal.h.

47 { return m_peak_distance; }

References m_peak_distance.

◆ positionVariance()

double IInterference::positionVariance ( ) const
inlineinherited

Returns the position variance.

Definition at line 40 of file IInterference.h.

40 { return m_position_var; }

References IInterference::m_position_var.

◆ randomSample()

double InterferenceRadialParaCrystal::randomSample ( ) const
inline

Definition at line 53 of file InterferenceRadialParaCrystal.h.

53 { return m_pdf->createSampler()->randomSample(); }

References m_pdf.

◆ setDomainSize()

void InterferenceRadialParaCrystal::setDomainSize ( double  size)

Sets domain size (finite size corrections).

Parameters
sizesize of coherence domain along the lattice main axis in nanometers

Definition at line 65 of file InterferenceRadialParaCrystal.cpp.

66 {
67  m_domain_size = size;
68 }

References m_domain_size.

◆ setKappa()

void InterferenceRadialParaCrystal::setKappa ( double  kappa)

Sets size spacing coupling parameter of the Size Spacing Correlation Approximation.

Definition at line 52 of file InterferenceRadialParaCrystal.cpp.

53 {
54  m_kappa = kappa;
55 }

References kappa(), and m_kappa.

Referenced by ExemplarySamples::createSizeDistributionSSCAModel().

Here is the call graph for this function:

◆ setPositionVariance()

void IInterference::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 IInterference.cpp.

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

References IInterference::m_position_var.

Referenced by ExemplarySamples::createFiniteSquareLattice2D(), and ExemplarySamples::createSuperLattice().

◆ setProbabilityDistribution()

void InterferenceRadialParaCrystal::setProbabilityDistribution ( const IProfile1D pdf)

Sets one-dimensional probability distribution.

Parameters
pdfprobability distribution (Fourier transform of probability density)

Definition at line 83 of file InterferenceRadialParaCrystal.cpp.

84 {
85  m_pdf.reset(pdf.clone());
86 }
IProfile1D * clone() const override=0

References IProfile1D::clone(), and m_pdf.

Referenced by ExemplarySamples::createCosineRipple(), ExemplarySamples::createRadialParaCrystal(), ExemplarySamples::createSizeDistributionDAModel(), ExemplarySamples::createSizeDistributionLMAModel(), ExemplarySamples::createSizeDistributionSSCAModel(), and ExemplarySamples::createTriangularRipple().

Here is the call graph for this function:

◆ structureFactor()

double IInterference::structureFactor ( R3  q,
double  outer_iff = 1.0 
) const
virtualinherited

The interference function for a given wavevector transfer.

Reimplemented in Interference2DSuperLattice.

Definition at line 35 of file IInterference.cpp.

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

References IInterference::iff_no_inner().

Here is the call graph for this function:

◆ supportsMultilayer()

virtual bool IInterference::supportsMultilayer ( ) const
inlinevirtualinherited

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

Reimplemented in InterferenceFinite3DLattice, and Interference3DLattice.

Definition at line 47 of file IInterference.h.

47 { return true; }

Referenced by IInterference::DWfactor().

◆ transferToCPP()

virtual void ICloneable::transferToCPP ( )
inlinevirtualinherited

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

Definition at line 32 of file ICloneable.h.

Member Data Documentation

◆ m_damping_length

double InterferenceRadialParaCrystal::m_damping_length
private

damping length of paracrystal Fourier transformed probability distribution of the nearest particle

Definition at line 60 of file InterferenceRadialParaCrystal.h.

Referenced by InterferenceRadialParaCrystal(), clone(), dampingLength(), and FTPDF().

◆ m_domain_size

double InterferenceRadialParaCrystal::m_domain_size
private

Size of coherence domain.

Definition at line 65 of file InterferenceRadialParaCrystal.h.

Referenced by InterferenceRadialParaCrystal(), clone(), domainSize(), iff_without_dw(), and setDomainSize().

◆ m_kappa

double InterferenceRadialParaCrystal::m_kappa
private

Size-spacing coupling parameter.

Definition at line 64 of file InterferenceRadialParaCrystal.h.

Referenced by InterferenceRadialParaCrystal(), clone(), kappa(), and setKappa().

◆ m_P

◆ m_pdf

std::unique_ptr<IProfile1D> InterferenceRadialParaCrystal::m_pdf
private

◆ m_peak_distance

double InterferenceRadialParaCrystal::m_peak_distance
private

the distance to the first neighbor peak

Definition at line 59 of file InterferenceRadialParaCrystal.h.

Referenced by InterferenceRadialParaCrystal(), clone(), FTPDF(), iff_without_dw(), and peakDistance().

◆ m_position_var

◆ m_use_damping_length

bool InterferenceRadialParaCrystal::m_use_damping_length
private

Definition at line 63 of file InterferenceRadialParaCrystal.h.

Referenced by InterferenceRadialParaCrystal(), and FTPDF().


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