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

Description

Interference function of a 2D lattice.

Definition at line 25 of file Interference2DLattice.h.

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

Public Member Functions

 Interference2DLattice (const Lattice2D &lattice)
 
 ~Interference2DLattice () override
 
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...
 
Interference2DLatticeclone () const override
 
double DWfactor (R3 q) const
 structureFactors the Debye-Waller factor for a given wavevector transfer More...
 
bool integrationOverXi () const
 
const Lattice2Dlattice () const
 
std::vector< const INode * > nodeChildren () const override
 Returns all children. More...
 
std::vector< const INode * > nodeOffspring () const
 Returns all descendants. More...
 
virtual std::vector< ParaMetaparDefs () const
 Returns the parameter definitions, to be hard-coded in each leaf class. More...
 
double particleDensity () const override
 Returns the particle density associated with this 2d lattice. More...
 
double positionVariance () const
 Returns the position variance. More...
 
void setDecayFunction (const IProfile2D &decay)
 Sets two-dimensional decay function. More...
 
void setIntegrationOverXi (bool integrate_xi)
 
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...
 
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

std::pair< double, double > calculateReciprocalVectorFraction (double qx, double qy, double xi) const
 Returns qx,qy coordinates of q - qint, where qint is a reciprocal lattice vector bounding the reciprocal unit cell to which q belongs. More...
 
double iff_without_dw (R3 q) const override
 Calculates the structure factor without Debye-Waller factor. More...
 
double interferenceAtOneRecLatticePoint (double qx, double qy) const
 Returns interference from a single reciprocal lattice vector. More...
 
double interferenceForXi (double xi, double qx, double qy) const
 
std::pair< double, double > rotateOrthonormal (double qx, double qy, double gamma) const
 Returns reciprocal coordinates in the coordinate system rotated by the angle gamma. More...
 

Private Attributes

std::unique_ptr< IProfile2Dm_decay
 
bool m_integrate_xi
 Integrate over the orientation xi. More...
 
std::unique_ptr< Lattice2Dm_lattice
 
int m_na
 
int m_nb
 determines the number of reciprocal lattice points to use More...
 
Lattice2D::ReciprocalBases m_sbase
 reciprocal lattice is stored without xi More...
 

Constructor & Destructor Documentation

◆ Interference2DLattice()

Interference2DLattice::Interference2DLattice ( const Lattice2D lattice)

Definition at line 50 of file Interference2DLattice.cpp.

51  : IInterference(0)
52  , m_integrate_xi(false)
53 {
54  m_lattice.reset(lattice.clone());
55  if (!m_lattice)
56  throw std::runtime_error("Interference2DLattice::initialize_rec_vectors() -> "
57  "Error. No lattice defined yet");
58 
59  BasicLattice2D base_lattice(m_lattice->length1(), m_lattice->length2(),
60  m_lattice->latticeAngle(), 0.);
61  m_sbase = base_lattice.reciprocalBases();
62 }
A two-dimensional Bravais lattice with no special symmetry.
Definition: Lattice2D.h:51
IInterference(const std::vector< double > &PValues)
std::unique_ptr< Lattice2D > m_lattice
bool m_integrate_xi
Integrate over the orientation xi.
const Lattice2D & lattice() const
Lattice2D::ReciprocalBases m_sbase
reciprocal lattice is stored without xi
Lattice2D * clone() const override=0

References Lattice2D::clone(), lattice(), m_lattice, m_sbase, and Lattice2D::reciprocalBases().

Referenced by clone().

Here is the call graph for this function:

◆ ~Interference2DLattice()

Interference2DLattice::~Interference2DLattice ( )
overridedefault

Member Function Documentation

◆ calculateReciprocalVectorFraction()

std::pair< double, double > Interference2DLattice::calculateReciprocalVectorFraction ( double  qx,
double  qy,
double  xi 
) const
private

Returns qx,qy coordinates of q - qint, where qint is a reciprocal lattice vector bounding the reciprocal unit cell to which q belongs.

Definition at line 171 of file Interference2DLattice.cpp.

172 {
173  double a = m_lattice->length1();
174  double b = m_lattice->length2();
175  double alpha = m_lattice->latticeAngle();
176  // first rotate the input to the system of m_sbase:
177  double qx_rot = qx * std::cos(xi) + qy * std::sin(xi);
178  double qy_rot = -qx * std::sin(xi) + qy * std::cos(xi);
179 
180  // find the reciprocal lattice coordinates of (qx_rot, qy_rot):
181  int qa_int = static_cast<int>(std::lround(a * qx_rot / M_TWOPI));
182  int qb_int = static_cast<int>(
183  std::lround(b * (qx_rot * std::cos(alpha) + qy_rot * std::sin(alpha)) / M_TWOPI));
184  // take the fractional part only (in m_sbase coordinates)
185  double qx_frac = qx_rot - qa_int * m_sbase.m_asx - qb_int * m_sbase.m_bsx;
186  double qy_frac = qy_rot - qa_int * m_sbase.m_asy - qb_int * m_sbase.m_bsy;
187  return {qx_frac, qy_frac};
188 }
#define M_TWOPI
Definition: Constants.h:54
double m_asy
x,y coordinates of a*
Definition: Lattice2D.h:30
double m_bsy
x,y coordinates of b*
Definition: Lattice2D.h:31

References Lattice2D::ReciprocalBases::m_asx, Lattice2D::ReciprocalBases::m_asy, Lattice2D::ReciprocalBases::m_bsx, Lattice2D::ReciprocalBases::m_bsy, m_lattice, m_sbase, and M_TWOPI.

Referenced by interferenceForXi().

◆ 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
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
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 Interference2DLattice::className ( ) const
inlinefinalvirtual

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

Implements INode.

Definition at line 31 of file Interference2DLattice.h.

31 { return "Interference2DLattice"; }

◆ clone()

Interference2DLattice * Interference2DLattice::clone ( ) const
overridevirtual

Implements IInterference.

Definition at line 66 of file Interference2DLattice.cpp.

67 {
68  auto* result = new Interference2DLattice(*m_lattice);
69  result->setPositionVariance(m_position_var);
70  result->setIntegrationOverXi(integrationOverXi());
71  if (m_decay)
72  result->setDecayFunction(*m_decay);
73  return result;
74 }
double m_position_var
Definition: IInterference.h:53
Interference2DLattice(const Lattice2D &lattice)
std::unique_ptr< IProfile2D > m_decay

References Interference2DLattice(), integrationOverXi(), m_decay, m_lattice, and IInterference::m_position_var.

Here is the call graph for this function:

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

◆ 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 Interference2DLattice::iff_without_dw ( R3  q) const
overrideprivatevirtual

Calculates the structure factor without Debye-Waller factor.

Implements IInterference.

Definition at line 120 of file Interference2DLattice.cpp.

121 {
122  if (!m_decay)
123  throw std::runtime_error("Interference2DLattice::evaluate"
124  " -> Error! No decay function defined.");
125  if (!m_integrate_xi)
126  return interferenceForXi(m_lattice->rotationAngle(), q.x(), q.y());
127  return RealIntegrator().integrate(
128  [&](double xi) -> double { return interferenceForXi(xi, q.x(), q.y()); }, 0.0,
129  M_TWOPI)
130  / M_TWOPI;
131 }
double interferenceForXi(double xi, double qx, double qy) const
To integrate a real function of a real variable.
Definition: IntegratorGK.h:28
double integrate(const std::function< double(double)> &f, double lmin, double lmax)

References RealIntegrator::integrate(), interferenceForXi(), m_decay, m_integrate_xi, m_lattice, and M_TWOPI.

Here is the call graph for this function:

◆ integrationOverXi()

bool Interference2DLattice::integrationOverXi ( ) const
inline

Definition at line 36 of file Interference2DLattice.h.

36 { return m_integrate_xi; }

References m_integrate_xi.

Referenced by clone().

◆ interferenceAtOneRecLatticePoint()

double Interference2DLattice::interferenceAtOneRecLatticePoint ( double  qx,
double  qy 
) const
private

Returns interference from a single reciprocal lattice vector.

Definition at line 148 of file Interference2DLattice.cpp.

149 {
150  if (!m_decay)
151  throw std::runtime_error("Interference2DLattice::interferenceAtOneRecLatticePoint"
152  " -> Error! No decay function defined.");
153  double gamma = m_decay->gamma();
154  auto qXY = rotateOrthonormal(qx, qy, gamma);
155  return m_decay->decayFT2D(qXY.first, qXY.second);
156 }
std::pair< double, double > rotateOrthonormal(double qx, double qy, double gamma) const
Returns reciprocal coordinates in the coordinate system rotated by the angle gamma.
double gamma(double x)

References ROOT::Math::Cephes::gamma(), m_decay, and rotateOrthonormal().

Referenced by interferenceForXi().

Here is the call graph for this function:

◆ interferenceForXi()

double Interference2DLattice::interferenceForXi ( double  xi,
double  qx,
double  qy 
) const
private

Definition at line 133 of file Interference2DLattice.cpp.

134 {
135  double result = 0.0;
136  auto q_frac = calculateReciprocalVectorFraction(qx, qy, xi);
137 
138  for (int i = -m_na - 1; i < m_na + 2; ++i) {
139  for (int j = -m_nb - 1; j < m_nb + 2; ++j) {
140  const double px = q_frac.first + i * m_sbase.m_asx + j * m_sbase.m_bsx;
141  const double py = q_frac.second + i * m_sbase.m_asy + j * m_sbase.m_bsy;
142  result += interferenceAtOneRecLatticePoint(px, py);
143  }
144  }
145  return particleDensity() * result;
146 }
std::pair< double, double > calculateReciprocalVectorFraction(double qx, double qy, double xi) const
Returns qx,qy coordinates of q - qint, where qint is a reciprocal lattice vector bounding the recipro...
int m_nb
determines the number of reciprocal lattice points to use
double interferenceAtOneRecLatticePoint(double qx, double qy) const
Returns interference from a single reciprocal lattice vector.
double particleDensity() const override
Returns the particle density associated with this 2d lattice.

References calculateReciprocalVectorFraction(), interferenceAtOneRecLatticePoint(), Lattice2D::ReciprocalBases::m_asx, Lattice2D::ReciprocalBases::m_asy, Lattice2D::ReciprocalBases::m_bsx, Lattice2D::ReciprocalBases::m_bsy, m_na, m_nb, m_sbase, and particleDensity().

Referenced by iff_without_dw().

Here is the call graph for this function:

◆ lattice()

const Lattice2D & Interference2DLattice::lattice ( ) const

Definition at line 101 of file Interference2DLattice.cpp.

102 {
103  if (!m_lattice)
104  throw std::runtime_error("Interference2DLattice::lattice() -> Error. "
105  "No lattice defined.");
106  return *m_lattice;
107 }

References m_lattice.

Referenced by Interference2DLattice().

◆ nodeChildren()

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

Returns all children.

Reimplemented from INode.

Definition at line 115 of file Interference2DLattice.cpp.

116 {
117  return std::vector<const INode*>() << m_decay << m_lattice;
118 }

References m_decay, and m_lattice.

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

virtual std::vector<ParaMeta> INode::parDefs ( ) const
inlinevirtualinherited

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

Reimplemented in ConstantBackground, GaussSphere, FuzzySphere, RotationEuler, RotationZ, RotationY, RotationX, Crystal, Layer, HexagonalLattice2D, SquareLattice2D, BasicLattice2D, LayerRoughness, TruncatedSpheroid, TruncatedSphere, TruncatedCube, Spheroid, Sphere, SawtoothRippleLorentz, SawtoothRippleGauss, SawtoothRippleBox, Pyramid6, Pyramid4, Pyramid3, Pyramid2, Prism6, Prism3, PlatonicTetrahedron, PlatonicOctahedron, LongBoxLorentz, LongBoxGauss, Icosahedron, HorizontalCylinder, HollowSphere, HemiEllipsoid, EllipsoidalCylinder, Dodecahedron, Cylinder, CosineRippleLorentz, CosineRippleGauss, CosineRippleBox, Cone, CantellatedCube, Box, Bipyramid4, BarLorentz, BarGauss, Profile2DVoigt, Profile2DCone, Profile2DGate, Profile2DGauss, Profile2DCauchy, Profile1DVoigt, Profile1DCosine, Profile1DTriangle, Profile1DGate, Profile1DGauss, Profile1DCauchy, MisesGaussPeakShape, MisesFisherGaussPeakShape, LorentzFisherPeakShape, GaussFisherPeakShape, IsotropicLorentzPeakShape, IsotropicGaussPeakShape, ParticleLayout, InterferenceTwin, InterferenceRadialParaCrystal, InterferenceHardDisk, Interference2DSuperLattice, Interference2DParaCrystal, Interference1DLattice, DistributionTrapezoid, DistributionCosine, DistributionLogNormal, DistributionGaussian, DistributionLorentz, DistributionGate, ResolutionFunction2DGaussian, PolFilter, FootprintSquare, and FootprintGauss.

Definition at line 51 of file INode.h.

51 { return {}; }

Referenced by INode::checkNodeArgs(), and IFormFactor::pythonConstructor().

◆ particleDensity()

double Interference2DLattice::particleDensity ( ) const
overridevirtual

Returns the particle density associated with this 2d lattice.

Reimplemented from IInterference.

Definition at line 109 of file Interference2DLattice.cpp.

110 {
111  double area = m_lattice->unitCellArea();
112  return area == 0.0 ? 0.0 : 1.0 / area;
113 }

References m_lattice.

Referenced by interferenceForXi().

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

◆ rotateOrthonormal()

std::pair< double, double > Interference2DLattice::rotateOrthonormal ( double  qx,
double  qy,
double  gamma 
) const
private

Returns reciprocal coordinates in the coordinate system rotated by the angle gamma.

Definition at line 159 of file Interference2DLattice.cpp.

161 {
162  double q_X = qx * std::cos(gamma) + qy * std::sin(gamma);
163  double q_Y = -qx * std::sin(gamma) + qy * std::cos(gamma);
164  return {q_X, q_Y};
165 }

References ROOT::Math::Cephes::gamma().

Referenced by interferenceAtOneRecLatticePoint().

Here is the call graph for this function:

◆ setDecayFunction()

void Interference2DLattice::setDecayFunction ( const IProfile2D decay)

Sets two-dimensional decay function.

Parameters
decaytwo-dimensional decay function in reciprocal space

Definition at line 78 of file Interference2DLattice.cpp.

79 {
80  m_decay.reset(decay.clone());
81  if (!m_decay)
82  throw std::runtime_error("Interference2DLattice::initialize_calc_factors"
83  " -> Error! No decay function defined.");
84 
85  // number of reciprocal lattice points to use
86  auto q_bounds = boundingReciprocalLatticeCoordinates(
87  nmax / m_decay->decayLengthX(), nmax / m_decay->decayLengthY(), m_lattice->length1(),
88  m_lattice->length2(), m_lattice->latticeAngle(), m_decay->gamma());
89  m_na = static_cast<int>(std::lround(q_bounds.first + 0.5));
90  m_nb = static_cast<int>(std::lround(q_bounds.second + 0.5));
91  m_na = std::max(m_na, min_points);
92  m_nb = std::max(m_nb, min_points);
93 }
IProfile2D * clone() const override=0

References IProfile2D::clone(), m_decay, m_lattice, m_na, and m_nb.

Referenced by ExemplarySamples::createBasic2DLattice(), ExemplarySamples::createBoxesSquareLattice2D(), ExemplarySamples::createCenteredSquareLattice2D(), ExemplarySamples::createParticleComposition(), ExemplarySamples::createRotatedSquareLattice2D(), and ExemplarySamples::createSquareLattice2D().

Here is the call graph for this function:

◆ setIntegrationOverXi()

void Interference2DLattice::setIntegrationOverXi ( bool  integrate_xi)

Definition at line 95 of file Interference2DLattice.cpp.

96 {
97  m_integrate_xi = integrate_xi;
98  m_lattice->setRotationEnabled(!m_integrate_xi); // deregister Xi in the case of integration
99 }

References m_integrate_xi, and m_lattice.

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

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

std::unique_ptr<IProfile2D> Interference2DLattice::m_decay
private

◆ m_integrate_xi

bool Interference2DLattice::m_integrate_xi
private

Integrate over the orientation xi.

Definition at line 61 of file Interference2DLattice.h.

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

◆ m_lattice

◆ m_na

int Interference2DLattice::m_na
private

Definition at line 65 of file Interference2DLattice.h.

Referenced by interferenceForXi(), and setDecayFunction().

◆ m_nb

int Interference2DLattice::m_nb
private

determines the number of reciprocal lattice points to use

Definition at line 65 of file Interference2DLattice.h.

Referenced by interferenceForXi(), and setDecayFunction().

◆ m_P

◆ m_position_var

◆ m_sbase

Lattice2D::ReciprocalBases Interference2DLattice::m_sbase
private

reciprocal lattice is stored without xi

Definition at line 64 of file Interference2DLattice.h.

Referenced by Interference2DLattice(), calculateReciprocalVectorFraction(), and interferenceForXi().


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