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

Public Member Functions

 Lattice ()
 
 Lattice (const kvector_t a1, const kvector_t a2, const kvector_t a3)
 
 Lattice (const Lattice &lattice)
 
 ~Lattice () override
 
void accept (INodeVisitor *visitor) const override
 
Lattice createTransformedLattice (const Transform3D &transform) const
 
void initialize () const
 
kvector_t getBasisVectorA () const
 
kvector_t getBasisVectorB () const
 
kvector_t getBasisVectorC () const
 
void resetBasis (const kvector_t a1, const kvector_t a2, const kvector_t a3)
 
kvector_t getMillerDirection (double h, double k, double l) const
 
double volume () const
 
void getReciprocalLatticeBasis (kvector_t &b1, kvector_t &b2, kvector_t &b3) const
 
ivector_t getNearestLatticeVectorCoordinates (const kvector_t vector_in) const
 
ivector_t getNearestReciprocalLatticeVectorCoordinates (const kvector_t vector_in) const
 
std::vector< kvector_treciprocalLatticeVectorsWithinRadius (const kvector_t input_vector, double radius) const
 
void setSelectionRule (const ISelectionRule &p_selection_rule)
 
void onChange () override
 
virtual std::string treeToString () const
 
void registerChild (INode *node)
 
virtual std::vector< const INode * > getChildren () const
 
virtual void setParent (const INode *newParent)
 
const INodeparent () const
 
INodeparent ()
 
int copyNumber (const INode *node) const
 
std::string displayName () const
 
ParameterPoolcreateParameterTree () const
 
ParameterPoolparameterPool () const
 
std::string parametersToString () const
 
RealParameterregisterParameter (const std::string &name, double *parpointer)
 
void registerVector (const std::string &base_name, kvector_t *p_vec, const std::string &units="nm")
 
void setParameterValue (const std::string &name, double value)
 
void setVectorValue (const std::string &base_name, kvector_t value)
 
RealParameterparameter (const std::string &name) const
 
void removeParameter (const std::string &name)
 
void removeVector (const std::string &base_name)
 
void setName (const std::string &name)
 
const std::string & getName () const
 

Static Public Member Functions

static Lattice createCubicLattice (double a)
 
static Lattice createFCCLattice (double a)
 
static Lattice createHexagonalLattice (double a, double c)
 
static Lattice createHCPLattice (double a, double c)
 
static Lattice createTetragonalLattice (double a, double c)
 
static Lattice createBCTLattice (double a, double c)
 
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 Attributes

const size_t m_NP
 
std::vector< double > m_P
 

Private Member Functions

Latticeoperator= (const Lattice &lattice)
 
void registerBasisVectors ()
 
std::vector< kvector_tvectorsWithinRadius (const kvector_t input_vector, const ivector_t &nearest_coords, double radius, const kvector_t v1, const kvector_t v2, const kvector_t v3, const kvector_t rec1, const kvector_t rec2, const kvector_t rec3) const
 
void computeReciprocalVectors () const
 

Static Private Member Functions

static void computeInverseVectors (const kvector_t v1, const kvector_t v2, const kvector_t v3, kvector_t o1, kvector_t o2, kvector_t o3)
 

Private Attributes

ISelectionRulemp_selection_rule
 
kvector_t m_a
 
kvector_t m_b
 
kvector_t m_c
 
kvector_t m_ra
 
kvector_t m_rb
 
kvector_t m_rc
 
bool m_cache_ok
 
const INodem_parent {nullptr}
 
std::string m_name
 
std::unique_ptr< ParameterPoolm_pool
 

Detailed Description

A lattice with three basis vectors.

Definition at line 27 of file Lattice.h.

Constructor & Destructor Documentation

◆ Lattice() [1/3]

Lattice::Lattice ( )

Definition at line 22 of file Lattice.cpp.

23  : mp_selection_rule(nullptr), m_a({1.0, 0.0, 0.0}), m_b({0.0, 1.0, 0.0}), m_c({0.0, 0.0, 1.0}),
24  m_cache_ok(false)
25 {
26  setName("Lattice");
27  initialize();
29 }
void setName(const std::string &name)
bool m_cache_ok
Boolean indicating if the reciprocal vectors are already initialized in the cache.
Definition: Lattice.h:116
void registerBasisVectors()
Definition: Lattice.cpp:190
kvector_t m_b
Definition: Lattice.h:113
kvector_t m_c
Basis vectors in real space.
Definition: Lattice.h:113
void initialize() const
Initializes cached data.
Definition: Lattice.cpp:66
ISelectionRule * mp_selection_rule
Definition: Lattice.h:112
kvector_t m_a
Definition: Lattice.h:113

Referenced by createBCTLattice(), createCubicLattice(), createFCCLattice(), createHCPLattice(), createHexagonalLattice(), and createTetragonalLattice().

◆ Lattice() [2/3]

Lattice::Lattice ( const kvector_t  a1,
const kvector_t  a2,
const kvector_t  a3 
)

Definition at line 31 of file Lattice.cpp.

32  : mp_selection_rule(nullptr), m_a(a1), m_b(a2), m_c(a3), m_cache_ok(false)
33 {
34  setName("Lattice");
35  initialize();
37 }

References initialize(), registerBasisVectors(), and IParameterized::setName().

Here is the call graph for this function:

◆ Lattice() [3/3]

Lattice::Lattice ( const Lattice lattice)

Definition at line 39 of file Lattice.cpp.

40  : INode(), mp_selection_rule(nullptr), m_a(lattice.m_a), m_b(lattice.m_b), m_c(lattice.m_c),
41  m_cache_ok(false)
42 {
43  setName("Lattice");
44  initialize();
45  if (lattice.mp_selection_rule)
48 }
INode()
Definition: INode.h:51
void setSelectionRule(const ISelectionRule &p_selection_rule)
Sets a selection rule for the reciprocal vectors.
Definition: Lattice.cpp:277

References initialize(), mp_selection_rule, registerBasisVectors(), IParameterized::setName(), and setSelectionRule().

Here is the call graph for this function:

◆ ~Lattice()

Lattice::~Lattice ( )
override

Definition at line 50 of file Lattice.cpp.

51 {
52  delete mp_selection_rule;
53 }

References mp_selection_rule.

Member Function Documentation

◆ accept()

void Lattice::accept ( INodeVisitor visitor) const
inlineoverridevirtual

Calls the INodeVisitor's visit method.

Implements INode.

Definition at line 35 of file Lattice.h.

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

References INodeVisitor::visit().

Here is the call graph for this function:

◆ createTransformedLattice()

Lattice Lattice::createTransformedLattice ( const Transform3D transform) const

Creates transformed lattice.

Definition at line 55 of file Lattice.cpp.

56 {
57  kvector_t a1 = transform.transformed(m_a);
58  kvector_t a2 = transform.transformed(m_b);
59  kvector_t a3 = transform.transformed(m_c);
60  Lattice result = {a1, a2, a3};
63  return result;
64 }
A lattice with three basis vectors.
Definition: Lattice.h:28
ivector_t transformed(const ivector_t &v) const
Return transformed vector v.

References m_a, m_b, m_c, mp_selection_rule, setSelectionRule(), and Transform3D::transformed().

Referenced by LatticeUtils::createBCTLattice(), LatticeUtils::createFCCLattice(), LatticeUtils::createHCPLattice(), and Crystal::transformedLattice().

Here is the call graph for this function:

◆ initialize()

void Lattice::initialize ( ) const

Initializes cached data.

Definition at line 66 of file Lattice.cpp.

67 {
69  m_cache_ok = true;
70 }
void computeReciprocalVectors() const
Definition: Lattice.cpp:199

References computeReciprocalVectors(), and m_cache_ok.

Referenced by getReciprocalLatticeBasis(), Lattice(), and reciprocalLatticeVectorsWithinRadius().

Here is the call graph for this function:

◆ getBasisVectorA()

kvector_t Lattice::getBasisVectorA ( ) const
inline

◆ getBasisVectorB()

kvector_t Lattice::getBasisVectorB ( ) const
inline

◆ getBasisVectorC()

kvector_t Lattice::getBasisVectorC ( ) const
inline

◆ resetBasis()

void Lattice::resetBasis ( const kvector_t  a1,
const kvector_t  a2,
const kvector_t  a3 
)

Resets the basis vectors.

Definition at line 72 of file Lattice.cpp.

73 {
74  m_a = a1;
75  m_b = a2;
76  m_c = a3;
77  onChange();
78 }
void onChange() override
Action to be taken in inherited class when a parameter has changed.
Definition: Lattice.cpp:185

References m_a, m_b, m_c, and onChange().

Referenced by MillerIndexOrientation::usePrimitiveLattice().

Here is the call graph for this function:

◆ getMillerDirection()

kvector_t Lattice::getMillerDirection ( double  h,
double  k,
double  l 
) const

Returns normalized direction corresponding to the given Miller indices.

Definition at line 80 of file Lattice.cpp.

81 {
82  kvector_t b1, b2, b3;
83  getReciprocalLatticeBasis(b1, b2, b3);
84  kvector_t direction = h * b1 + k * b2 + l * b3;
85  return direction.unit();
86 }
BasicVector3D< T > unit() const
Returns unit vector in direction of this. Throws for null vector.
void getReciprocalLatticeBasis(kvector_t &b1, kvector_t &b2, kvector_t &b3) const
Returns the reciprocal basis vectors.
Definition: Lattice.cpp:93

References getReciprocalLatticeBasis(), and BasicVector3D< T >::unit().

Referenced by MillerIndexOrientation::transformationMatrix().

Here is the call graph for this function:

◆ volume()

double Lattice::volume ( ) const

Returns the volume of the unit cell.

Definition at line 88 of file Lattice.cpp.

89 {
90  return std::abs(m_a.dot(m_b.cross(m_c)));
91 }
auto dot(const BasicVector3D< U > &v) const
Returns dot product of vectors (antilinear in the first [=self] argument).
auto cross(const BasicVector3D< U > &v) const
Returns cross product of vectors (linear in both arguments).

References BasicVector3D< T >::cross(), BasicVector3D< T >::dot(), m_a, m_b, and m_c.

Referenced by FormFactorCrystal::evaluate(), FormFactorCrystal::evaluatePol(), and Crystal::homogeneousRegions().

Here is the call graph for this function:

◆ getReciprocalLatticeBasis()

void Lattice::getReciprocalLatticeBasis ( kvector_t b1,
kvector_t b2,
kvector_t b3 
) const

Returns the reciprocal basis vectors.

Definition at line 93 of file Lattice.cpp.

94 {
95  if (!m_cache_ok) {
96  initialize();
97  }
98  b1 = m_ra;
99  b2 = m_rb;
100  b3 = m_rc;
101  return;
102 }
kvector_t m_rc
Cache of basis vectors in reciprocal space.
Definition: Lattice.h:114
kvector_t m_ra
Definition: Lattice.h:114
kvector_t m_rb
Definition: Lattice.h:114

References initialize(), m_cache_ok, m_ra, m_rb, and m_rc.

Referenced by getMillerDirection().

Here is the call graph for this function:

◆ getNearestLatticeVectorCoordinates()

ivector_t Lattice::getNearestLatticeVectorCoordinates ( const kvector_t  vector_in) const

Returns the nearest lattice point from a given vector.

Definition at line 104 of file Lattice.cpp.

105 {
106  double a1_coord = vector_in.dot(m_ra) / M_TWOPI;
107  double a2_coord = vector_in.dot(m_rb) / M_TWOPI;
108  double a3_coord = vector_in.dot(m_rc) / M_TWOPI;
109  int c1 = static_cast<int>(std::floor(a1_coord + 0.5));
110  int c2 = static_cast<int>(std::floor(a2_coord + 0.5));
111  int c3 = static_cast<int>(std::floor(a3_coord + 0.5));
112  return ivector_t(c1, c2, c3);
113 }
#define M_TWOPI
Definition: MathConstants.h:49
BasicVector3D< int > ivector_t
Definition: Vectors3D.h:20

References BasicVector3D< T >::dot(), m_ra, m_rb, m_rc, and M_TWOPI.

Here is the call graph for this function:

◆ getNearestReciprocalLatticeVectorCoordinates()

ivector_t Lattice::getNearestReciprocalLatticeVectorCoordinates ( const kvector_t  vector_in) const

Returns the nearest reciprocal lattice point from a given vector.

Definition at line 115 of file Lattice.cpp.

116 {
117  double b1_coord = vector_in.dot(m_a) / M_TWOPI;
118  double b2_coord = vector_in.dot(m_b) / M_TWOPI;
119  double b3_coord = vector_in.dot(m_c) / M_TWOPI;
120  int c1 = static_cast<int>(std::floor(b1_coord + 0.5));
121  int c2 = static_cast<int>(std::floor(b2_coord + 0.5));
122  int c3 = static_cast<int>(std::floor(b3_coord + 0.5));
123  return ivector_t(c1, c2, c3);
124 }

References BasicVector3D< T >::dot(), m_a, m_b, m_c, and M_TWOPI.

Referenced by reciprocalLatticeVectorsWithinRadius().

Here is the call graph for this function:

◆ reciprocalLatticeVectorsWithinRadius()

std::vector< kvector_t > Lattice::reciprocalLatticeVectorsWithinRadius ( const kvector_t  input_vector,
double  radius 
) const

Computes a list of reciprocal lattice vectors within a specified distance of a given vector.

Definition at line 126 of file Lattice.cpp.

128 {
129  if (!m_cache_ok)
130  initialize();
131  ivector_t nearest_coords = getNearestReciprocalLatticeVectorCoordinates(input_vector);
132  return vectorsWithinRadius(input_vector, nearest_coords, radius, m_ra, m_rb, m_rc, m_a, m_b,
133  m_c);
134 }
std::vector< kvector_t > vectorsWithinRadius(const kvector_t input_vector, const ivector_t &nearest_coords, double radius, const kvector_t v1, const kvector_t v2, const kvector_t v3, const kvector_t rec1, const kvector_t rec2, const kvector_t rec3) const
Definition: Lattice.cpp:209
ivector_t getNearestReciprocalLatticeVectorCoordinates(const kvector_t vector_in) const
Returns the nearest reciprocal lattice point from a given vector.
Definition: Lattice.cpp:115
const double radius(5 *Units::nanometer)

References getNearestReciprocalLatticeVectorCoordinates(), initialize(), m_a, m_b, m_c, m_cache_ok, m_ra, m_rb, m_rc, anonymous_namespace{SlicedCylindersBuilder.cpp}::radius(), and vectorsWithinRadius().

Referenced by FormFactorCrystal::evaluate(), FormFactorCrystal::evaluatePol(), and InterferenceFunction3DLattice::iff_without_dw().

Here is the call graph for this function:

◆ setSelectionRule()

void Lattice::setSelectionRule ( const ISelectionRule p_selection_rule)

Sets a selection rule for the reciprocal vectors.

Definition at line 277 of file Lattice.cpp.

278 {
279  delete mp_selection_rule;
280  mp_selection_rule = p_selection_rule.clone();
281 }
virtual ISelectionRule * clone() const =0

References ISelectionRule::clone(), and mp_selection_rule.

Referenced by createTransformedLattice(), and Lattice().

Here is the call graph for this function:

◆ createCubicLattice()

Lattice Lattice::createCubicLattice ( double  a)
static

Returns a primitive cubic (cP) lattice with edge length a.

Definition at line 136 of file Lattice.cpp.

137 {
138  kvector_t a1(a, 0.0, 0.0);
139  kvector_t a2(0.0, a, 0.0);
140  kvector_t a3(0.0, 0.0, a);
141  return Lattice(a1, a2, a3);
142 }
Lattice()
Definition: Lattice.cpp:22

References Lattice().

Referenced by LatticeUtils::createFCCLattice().

Here is the call graph for this function:

◆ createFCCLattice()

Lattice Lattice::createFCCLattice ( double  a)
static

Returns a face-centered cubic (cF) lattice with edge length a.

Definition at line 144 of file Lattice.cpp.

145 {
146  double b = a / 2.0;
147  kvector_t a1(0.0, b, b);
148  kvector_t a2(b, 0.0, b);
149  kvector_t a3(b, b, 0.0);
150  return Lattice(a1, a2, a3);
151 }

References Lattice().

Referenced by LatticeUtils::createFCCLattice().

Here is the call graph for this function:

◆ createHexagonalLattice()

Lattice Lattice::createHexagonalLattice ( double  a,
double  c 
)
static

Returns a primitive hexagonal (hP) lattice with hexagonal edge a and height c.

Definition at line 153 of file Lattice.cpp.

154 {
155  kvector_t a1(a, 0.0, 0.0);
156  kvector_t a2(-a / 2.0, std::sqrt(3.0) * a / 2.0, 0.0);
157  kvector_t a3(0.0, 0.0, c);
158  return Lattice(a1, a2, a3);
159 }

References Lattice().

Referenced by LatticeUtils::createHCPLattice().

Here is the call graph for this function:

◆ createHCPLattice()

Lattice Lattice::createHCPLattice ( double  a,
double  c 
)
static

TODO: Clarify how this is meant: HCP is not a Bravais lattice.

Definition at line 161 of file Lattice.cpp.

162 {
163  kvector_t a1(a, 0.0, 0.0);
164  kvector_t a2(-a / 2.0, std::sqrt(3.0) * a / 2.0, 0);
165  kvector_t a3(a / 2.0, a / std::sqrt(3.0) / 2.0, c / 2.0);
166  return Lattice(a1, a2, a3);
167 }

References Lattice().

Referenced by LatticeUtils::createHCPLattice().

Here is the call graph for this function:

◆ createTetragonalLattice()

Lattice Lattice::createTetragonalLattice ( double  a,
double  c 
)
static

Returns a primitive tetragonal (tP) lattice with square base edge a and height c.

Definition at line 169 of file Lattice.cpp.

170 {
171  kvector_t a1(a, 0.0, 0.0);
172  kvector_t a2(0.0, a, 0.0);
173  kvector_t a3(0.0, 0.0, c);
174  return Lattice(a1, a2, a3);
175 }

References Lattice().

Referenced by LatticeUtils::createBCTLattice().

Here is the call graph for this function:

◆ createBCTLattice()

Lattice Lattice::createBCTLattice ( double  a,
double  c 
)
static

Returns a body-centered cubic (cI) lattice with edge length a.

TODO: Clarify meaning of c

Definition at line 177 of file Lattice.cpp.

178 {
179  kvector_t a1(a, 0.0, 0.0);
180  kvector_t a2(0.0, a, 0.0);
181  kvector_t a3(a / 2.0, a / 2.0, c / 2.0);
182  return Lattice(a1, a2, a3);
183 }

References Lattice().

Referenced by LatticeUtils::createBCTLattice().

Here is the call graph for this function:

◆ onChange()

void Lattice::onChange ( )
overridevirtual

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

Reimplemented from IParameterized.

Definition at line 185 of file Lattice.cpp.

186 {
187  m_cache_ok = false;
188 }

References m_cache_ok.

Referenced by resetBasis().

◆ operator=()

Lattice& Lattice::operator= ( const Lattice lattice)
private

◆ registerBasisVectors()

void Lattice::registerBasisVectors ( )
private

Definition at line 190 of file Lattice.cpp.

191 {
192  if (!parameter(XComponentName("BasisA"))) {
193  registerVector("BasisA", &m_a, "nm");
194  registerVector("BasisB", &m_b, "nm");
195  registerVector("BasisC", &m_c, "nm");
196  }
197 }
RealParameter * parameter(const std::string &name) const
Returns parameter with given 'name'.
static std::string XComponentName(const std::string &base_name)
void registerVector(const std::string &base_name, kvector_t *p_vec, const std::string &units="nm")

References m_a, m_b, m_c, IParameterized::parameter(), IParameterized::registerVector(), and IParameterized::XComponentName().

Referenced by Lattice().

Here is the call graph for this function:

◆ vectorsWithinRadius()

std::vector< kvector_t > Lattice::vectorsWithinRadius ( const kvector_t  input_vector,
const ivector_t nearest_coords,
double  radius,
const kvector_t  v1,
const kvector_t  v2,
const kvector_t  v3,
const kvector_t  rec1,
const kvector_t  rec2,
const kvector_t  rec3 
) const
private

Definition at line 209 of file Lattice.cpp.

215 {
216  int max_X = static_cast<int>(std::floor(rec1.mag() * radius / M_TWOPI + 0.5));
217  int max_Y = static_cast<int>(std::floor(rec2.mag() * radius / M_TWOPI + 0.5));
218  int max_Z = static_cast<int>(std::floor(rec3.mag() * radius / M_TWOPI + 0.5));
219 
220  std::vector<kvector_t> ret;
221  for (int index_X = -max_X; index_X <= max_X; ++index_X) {
222  for (int index_Y = -max_Y; index_Y <= max_Y; ++index_Y) {
223  for (int index_Z = -max_Z; index_Z <= max_Z; ++index_Z) {
224  ivector_t coords(index_X + nearest_coords[0], index_Y + nearest_coords[1],
225  index_Z + nearest_coords[2]);
227  continue;
228  kvector_t latticePoint = coords[0] * v1 + coords[1] * v2 + coords[2] * v3;
229  if ((latticePoint - input_vector).mag() <= radius)
230  ret.push_back(latticePoint);
231  }
232  }
233  }
234  return ret;
235 }
double mag() const
Returns magnitude of the vector.
virtual bool coordinateSelected(const ivector_t &coordinate) const =0

References ISelectionRule::coordinateSelected(), M_TWOPI, BasicVector3D< T >::mag(), mp_selection_rule, and anonymous_namespace{SlicedCylindersBuilder.cpp}::radius().

Referenced by reciprocalLatticeVectorsWithinRadius().

Here is the call graph for this function:

◆ computeReciprocalVectors()

void Lattice::computeReciprocalVectors ( ) const
private

Definition at line 199 of file Lattice.cpp.

200 {
201  kvector_t a23 = m_b.cross(m_c);
202  kvector_t a31 = m_c.cross(m_a);
203  kvector_t a12 = m_a.cross(m_b);
204  m_ra = M_TWOPI / m_a.dot(a23) * a23;
205  m_rb = M_TWOPI / m_b.dot(a31) * a31;
206  m_rc = M_TWOPI / m_c.dot(a12) * a12;
207 }

References BasicVector3D< T >::cross(), BasicVector3D< T >::dot(), m_a, m_b, m_c, m_ra, m_rb, m_rc, and M_TWOPI.

Referenced by initialize().

Here is the call graph for this function:

◆ computeInverseVectors()

void Lattice::computeInverseVectors ( const kvector_t  v1,
const kvector_t  v2,
const kvector_t  v3,
kvector_t  o1,
kvector_t  o2,
kvector_t  o3 
)
staticprivate

Definition at line 237 of file Lattice.cpp.

239 {
240  gsl_matrix* p_basisMatrix = gsl_matrix_alloc(3, 3);
241  gsl_matrix* p_inverseMatrix = gsl_matrix_alloc(3, 3);
242  gsl_permutation* p_perm = gsl_permutation_alloc(3);
243  int s;
244 
245  gsl_matrix_set(p_basisMatrix, 0, 0, v1.x());
246  gsl_matrix_set(p_basisMatrix, 0, 1, v2.x());
247  gsl_matrix_set(p_basisMatrix, 0, 2, v3.x());
248 
249  gsl_matrix_set(p_basisMatrix, 1, 0, v1.y());
250  gsl_matrix_set(p_basisMatrix, 1, 1, v2.y());
251  gsl_matrix_set(p_basisMatrix, 1, 2, v3.y());
252 
253  gsl_matrix_set(p_basisMatrix, 2, 0, v1.z());
254  gsl_matrix_set(p_basisMatrix, 2, 1, v2.z());
255  gsl_matrix_set(p_basisMatrix, 2, 2, v3.z());
256 
257  gsl_linalg_LU_decomp(p_basisMatrix, p_perm, &s);
258  gsl_linalg_LU_invert(p_basisMatrix, p_perm, p_inverseMatrix);
259 
260  o1.setX(gsl_matrix_get(p_inverseMatrix, 0, 0));
261  o1.setY(gsl_matrix_get(p_inverseMatrix, 1, 0));
262  o1.setZ(gsl_matrix_get(p_inverseMatrix, 2, 0));
263 
264  o2.setX(gsl_matrix_get(p_inverseMatrix, 0, 1));
265  o2.setY(gsl_matrix_get(p_inverseMatrix, 1, 1));
266  o2.setZ(gsl_matrix_get(p_inverseMatrix, 2, 1));
267 
268  o3.setX(gsl_matrix_get(p_inverseMatrix, 0, 2));
269  o3.setY(gsl_matrix_get(p_inverseMatrix, 1, 2));
270  o3.setZ(gsl_matrix_get(p_inverseMatrix, 2, 2));
271 
272  gsl_permutation_free(p_perm);
273  gsl_matrix_free(p_basisMatrix);
274  gsl_matrix_free(p_inverseMatrix);
275 }
void setX(const T &a)
Sets x-component in cartesian coordinate system.
Definition: BasicVector3D.h:71
T z() const
Returns z-component in cartesian coordinate system.
Definition: BasicVector3D.h:68
void setY(const T &a)
Sets y-component in cartesian coordinate system.
Definition: BasicVector3D.h:73
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
void setZ(const T &a)
Sets z-component in cartesian coordinate system.
Definition: BasicVector3D.h:75

References BasicVector3D< T >::setX(), BasicVector3D< T >::setY(), BasicVector3D< T >::setZ(), BasicVector3D< T >::x(), BasicVector3D< T >::y(), and BasicVector3D< T >::z().

Here is the call graph for this function:

◆ 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(), InterferenceFunction2DParaCrystal::setLattice(), InterferenceFunction2DSuperLattice::setLattice(), InterferenceFunctionFinite2DLattice::setLattice(), InterferenceFunctionRadialParaCrystal::setProbabilityDistribution(), InterferenceFunction2DParaCrystal::setProbabilityDistributions(), ConvolutionDetectorResolution::setResolutionFunction(), IParticle::setRotation(), LayerInterface::setRoughness(), and InterferenceFunction2DSuperLattice::setSubstructureIFF().

Here is the call graph for this function:

◆ getChildren()

◆ setParent()

void INode::setParent ( const INode newParent)
virtualinherited

Reimplemented in SampleProvider.

Definition at line 69 of file INode.cpp.

70 {
71  m_parent = newParent;
72 }
const INode * m_parent
Definition: INode.h:81

References INode::m_parent.

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

◆ parent() [1/2]

const INode * INode::parent ( ) const
inherited

◆ parent() [2/2]

INode * INode::parent ( )
inherited

Definition at line 79 of file INode.cpp.

80 {
81  return const_cast<INode*>(m_parent);
82 }
Base class for tree-like structures containing parameterized objects.
Definition: INode.h:49

References INode::m_parent.

◆ copyNumber()

int INode::copyNumber ( const INode node) const
inherited

Returns copyNumber of child, which takes into account existence of children with same name.

Definition at line 84 of file INode.cpp.

85 {
86  if (node->parent() != this)
87  return -1;
88 
89  int result(-1), count(0);
90  for (auto child : getChildren()) {
91 
92  if (child == nullptr)
93  throw std::runtime_error("INode::copyNumber() -> Error. Nullptr as child.");
94 
95  if (child == node)
96  result = count;
97 
98  if (child->getName() == node->getName())
99  ++count;
100  }
101 
102  return count > 1 ? result : -1;
103 }
const INode * parent() const
Definition: INode.cpp:74
virtual std::vector< const INode * > getChildren() const
Returns a vector of children (const).
Definition: INode.cpp:64
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::InterferenceFunction2DParaCrystal(), InterferenceFunctionHardDisk::InterferenceFunctionHardDisk(), InterferenceFunctionRadialParaCrystal::InterferenceFunctionRadialParaCrystal(), InterferenceFunctionTwin::InterferenceFunctionTwin(), Lattice2D::Lattice2D(), LayerRoughness::LayerRoughness(), MultiLayer::MultiLayer(), ParticleDistribution::ParticleDistribution(), PlainMultiLayerBySLDBuilder::PlainMultiLayerBySLDBuilder(), IParticle::registerAbundance(), ParticleLayout::registerParticleDensity(), Layer::registerThickness(), IParameterized::registerVector(), ParticleLayout::registerWeight(), ResolutionFunction2DGaussian::ResolutionFunction2DGaussian(), ResonatorBuilder::ResonatorBuilder(), Lattice2D::setRotationEnabled(), SquareLattice::SquareLattice(), and TriangularRippleBuilder::TriangularRippleBuilder().

Here is the call graph for this function:

◆ registerVector()

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

Definition at line 54 of file IParameterized.cpp.

56 {
57  registerParameter(XComponentName(base_name), &((*p_vec)[0])).setUnit(units);
58  registerParameter(YComponentName(base_name), &((*p_vec)[1])).setUnit(units);
59  registerParameter(ZComponentName(base_name), &((*p_vec)[2])).setUnit(units);
60 }
RealParameter & registerParameter(const std::string &name, double *parpointer)
static std::string YComponentName(const std::string &base_name)
static std::string ZComponentName(const std::string &base_name)
RealParameter & setUnit(const std::string &name)

References IParameterized::registerParameter(), RealParameter::setUnit(), IParameterized::XComponentName(), IParameterized::YComponentName(), and IParameterized::ZComponentName().

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

Here is the call graph for this function:

◆ setParameterValue()

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

Definition at line 62 of file IParameterized.cpp.

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

References IParameterized::createParameterTree(), IParameterized::m_pool, and ParameterPool::setMatchedParametersValue().

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

Here is the call graph for this function:

◆ setVectorValue()

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

Definition at line 75 of file IParameterized.cpp.

76 {
77  setParameterValue(XComponentName(base_name), value.x());
78  setParameterValue(YComponentName(base_name), value.y());
79  setParameterValue(ZComponentName(base_name), value.z());
80 }
void setParameterValue(const std::string &name, double value)

References IParameterized::setParameterValue(), BasicVector3D< T >::x(), IParameterized::XComponentName(), BasicVector3D< T >::y(), IParameterized::YComponentName(), BasicVector3D< T >::z(), and IParameterized::ZComponentName().

Here is the call graph for this function:

◆ parameter()

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

◆ removeParameter()

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

◆ removeVector()

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

Definition at line 93 of file IParameterized.cpp.

94 {
95  removeParameter(XComponentName(base_name));
96  removeParameter(YComponentName(base_name));
97  removeParameter(ZComponentName(base_name));
98 }
void removeParameter(const std::string &name)

References IParameterized::removeParameter(), IParameterized::XComponentName(), IParameterized::YComponentName(), and IParameterized::ZComponentName().

Referenced by IParticle::registerPosition().

Here is the call graph for this function:

◆ XComponentName()

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

◆ YComponentName()

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

Definition at line 105 of file IParameterized.cpp.

106 {
107  return base_name + "Y";
108 }

Referenced by IParameterized::registerVector(), IParameterized::removeVector(), and IParameterized::setVectorValue().

◆ ZComponentName()

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

Definition at line 110 of file IParameterized.cpp.

111 {
112  return base_name + "Z";
113 }

Referenced by IParameterized::registerVector(), IParameterized::removeVector(), and IParameterized::setVectorValue().

◆ setName()

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

Definition at line 68 of file IParameterized.h.

68 { m_name = name; }
std::string m_name

References IParameterized::m_name.

Referenced by BasicLattice::BasicLattice(), Beam::Beam(), Layer::clone(), ConvolutionDetectorResolution::ConvolutionDetectorResolution(), LayersWithAbsorptionBuilder::createSampleByIndex(), Basic2DParaCrystalBuilder::createSampleByIndex(), ParticleInVacuumBuilder::createSampleByIndex(), SimpleMagneticRotationBuilder::createSampleByIndex(), Crystal::Crystal(), DetectionProperties::DetectionProperties(), DistributionHandler::DistributionHandler(), FormFactorBAPol::FormFactorBAPol(), FormFactorCoreShell::FormFactorCoreShell(), FormFactorCrystal::FormFactorCrystal(), FormFactorDecoratorMaterial::FormFactorDecoratorMaterial(), FormFactorDecoratorPositionFactor::FormFactorDecoratorPositionFactor(), FormFactorDecoratorRotation::FormFactorDecoratorRotation(), FormFactorDWBA::FormFactorDWBA(), FormFactorDWBAPol::FormFactorDWBAPol(), FormFactorWeighted::FormFactorWeighted(), HexagonalLattice::HexagonalLattice(), IDetector::IDetector(), DepthProbeSimulation::initialize(), GISASSimulation::initialize(), OffSpecSimulation::initialize(), SpecularSimulation::initialize(), SpecularDetector1D::initialize(), MesoCrystal::initialize(), Particle::initialize(), ParticleComposition::initialize(), INode::INode(), Instrument::Instrument(), InterferenceFunction1DLattice::InterferenceFunction1DLattice(), InterferenceFunction2DLattice::InterferenceFunction2DLattice(), InterferenceFunction2DParaCrystal::InterferenceFunction2DParaCrystal(), InterferenceFunction2DSuperLattice::InterferenceFunction2DSuperLattice(), InterferenceFunction3DLattice::InterferenceFunction3DLattice(), InterferenceFunctionFinite2DLattice::InterferenceFunctionFinite2DLattice(), InterferenceFunctionFinite3DLattice::InterferenceFunctionFinite3DLattice(), InterferenceFunctionHardDisk::InterferenceFunctionHardDisk(), InterferenceFunctionNone::InterferenceFunctionNone(), InterferenceFunctionRadialParaCrystal::InterferenceFunctionRadialParaCrystal(), InterferenceFunctionTwin::InterferenceFunctionTwin(), ISampleBuilder::ISampleBuilder(), IsGISAXSDetector::IsGISAXSDetector(), Lattice(), 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

◆ mp_selection_rule

ISelectionRule* Lattice::mp_selection_rule
private

◆ m_a

◆ m_b

◆ m_c

◆ m_ra

◆ m_rb

◆ m_rc

kvector_t Lattice::m_rc
private

Cache of basis vectors in reciprocal space.

Definition at line 114 of file Lattice.h.

Referenced by computeReciprocalVectors(), getNearestLatticeVectorCoordinates(), getReciprocalLatticeBasis(), and reciprocalLatticeVectorsWithinRadius().

◆ m_cache_ok

bool Lattice::m_cache_ok
mutableprivate

Boolean indicating if the reciprocal vectors are already initialized in the cache.

Definition at line 116 of file Lattice.h.

Referenced by getReciprocalLatticeBasis(), initialize(), onChange(), and reciprocalLatticeVectorsWithinRadius().

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