BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
Lattice3D.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Lattice/Lattice3D.cpp
6 //! @brief Implements class Lattice.
7 //!
8 //! @homepage http://www.bornagainproject.org
9 //! @license GNU General Public License v3 or higher (see COPYING)
10 //! @copyright Forschungszentrum Jülich GmbH 2018
11 //! @authors Scientific Computing Group at MLZ (see CITATION, AUTHORS)
12 //
13 // ************************************************************************************************
14 
16 #include "Base/Math/Constants.h"
20 #include <gsl/gsl_linalg.h>
21 
23  : m_a(a), m_b(b), m_c(c)
24 {
25  setName("Lattice");
26  initialize();
27 }
28 
30  : INode(), m_a(lattice.m_a), m_b(lattice.m_b), m_c(lattice.m_c)
31 {
32  setName("Lattice");
33  initialize();
34  if (lattice.m_selection_rule)
36 }
37 
38 Lattice3D::~Lattice3D() = default;
39 
41 {
43  if (!parameter(XComponentName("BasisA"))) {
44  registerVector("BasisA", &m_a, "nm");
45  registerVector("BasisB", &m_b, "nm");
46  registerVector("BasisC", &m_c, "nm");
47  }
48 }
49 
51 {
53 }
54 
56 {
57  kvector_t q1 = transform.transformed(m_a);
58  kvector_t q2 = transform.transformed(m_b);
59  kvector_t q3 = transform.transformed(m_c);
60  Lattice3D result = {q1, q2, q3};
61  if (m_selection_rule)
63  return result;
64 }
65 
66 //! Currently unused but may be useful for checks
67 kvector_t Lattice3D::getMillerDirection(double h, double k, double l) const
68 {
69  kvector_t direction = h * m_ra + k * m_rb + l * m_rc;
70  return direction.unit();
71 }
72 
74 {
75  return std::abs(m_a.dot(m_b.cross(m_c)));
76 }
77 
78 //! Currently only used in tests
80 {
81  ra = m_ra;
82  rb = m_rb;
83  rc = m_rc;
84 }
85 
87 {
88  return {(int)std::lround(q.dot(m_a) / M_TWOPI), (int)std::lround(q.dot(m_b) / M_TWOPI),
89  (int)std::lround(q.dot(m_c) / M_TWOPI)};
90 }
91 
93  double dq) const
94 {
96 
97  int max_X = std::lround(m_a.mag() * dq / M_TWOPI);
98  int max_Y = std::lround(m_b.mag() * dq / M_TWOPI);
99  int max_Z = std::lround(m_c.mag() * dq / M_TWOPI);
100 
101  std::vector<kvector_t> ret;
102  for (int index_X = -max_X; index_X <= max_X; ++index_X) {
103  for (int index_Y = -max_Y; index_Y <= max_Y; ++index_Y) {
104  for (int index_Z = -max_Z; index_Z <= max_Z; ++index_Z) {
105  ivector_t coords(index_X + nearest_coords[0], index_Y + nearest_coords[1],
106  index_Z + nearest_coords[2]);
107  if (m_selection_rule && !m_selection_rule->coordinateSelected(coords))
108  continue;
109  kvector_t latticePoint = coords[0] * m_ra + coords[1] * m_rb + coords[2] * m_rc;
110  if ((latticePoint - q).mag() <= dq)
111  ret.push_back(latticePoint);
112  }
113  }
114  }
115  return ret;
116 }
117 
119 {
120  kvector_t q23 = m_b.cross(m_c);
121  kvector_t q31 = m_c.cross(m_a);
122  kvector_t q12 = m_a.cross(m_b);
123  m_ra = M_TWOPI / m_a.dot(q23) * q23;
124  m_rb = M_TWOPI / m_b.dot(q31) * q31;
125  m_rc = M_TWOPI / m_c.dot(q12) * q12;
126 }
127 
128 void Lattice3D::setSelectionRule(const ISelectionRule& selection_rule)
129 {
130  m_selection_rule.reset(selection_rule.clone());
131 }
Defines M_PI and some more mathematical constants.
#define M_TWOPI
Definition: Constants.h:54
Defines classes ISelectionRule, SimpleSelectionRule.
Defines class Lattice.
Defines class RealParameter.
Declares class Transform3D.
auto dot(const BasicVector3D< U > &v) const
Returns dot product of vectors (antilinear in the first [=self] argument).
BasicVector3D< T > unit() const
Returns unit vector in direction of this. Throws for null vector.
double mag() const
Returns magnitude of the vector.
auto cross(const BasicVector3D< U > &v) const
Returns cross product of vectors (linear in both arguments).
Base class for tree-like structures containing parameterized objects.
Definition: INode.h:49
void registerVector(const std::string &base_name, kvector_t *p_vec, const std::string &units="nm")
static std::string XComponentName(const std::string &base_name)
RealParameter * parameter(const std::string &name) const
Returns parameter with given 'name'.
void setName(const std::string &name)
Abstract base class for selection rules.
virtual ISelectionRule * clone() const =0
A Bravais lattice, characterized by three basis vectors, and optionally an ISelectionRule.
Definition: Lattice3D.h:29
kvector_t m_ra
Definition: Lattice3D.h:80
kvector_t m_rb
Definition: Lattice3D.h:80
double unitCellVolume() const
Returns the volume of the unit cell.
Definition: Lattice3D.cpp:73
ivector_t getNearestReciprocalLatticeVectorCoordinates(const kvector_t q) const
Returns the nearest reciprocal lattice point from a given vector.
Definition: Lattice3D.cpp:86
void computeReciprocalVectors() const
Definition: Lattice3D.cpp:118
kvector_t m_rc
Cache of basis vectors in reciprocal space.
Definition: Lattice3D.h:80
kvector_t m_c
Basis vectors in real space.
Definition: Lattice3D.h:77
void setSelectionRule(const ISelectionRule &selection_rule)
Sets a selection rule for the reciprocal vectors.
Definition: Lattice3D.cpp:128
Lattice3D transformed(const Transform3D &transform) const
Creates transformed lattice.
Definition: Lattice3D.cpp:55
void getReciprocalLatticeBasis(kvector_t &ra, kvector_t &rb, kvector_t &rc) const
Returns the reciprocal basis vectors.
Definition: Lattice3D.cpp:79
kvector_t getMillerDirection(double h, double k, double l) const
Returns normalized direction corresponding to the given Miller indices.
Definition: Lattice3D.cpp:67
std::unique_ptr< ISelectionRule > m_selection_rule
Definition: Lattice3D.h:78
kvector_t m_a
Definition: Lattice3D.h:77
std::vector< kvector_t > reciprocalLatticeVectorsWithinRadius(const kvector_t q, double dq) const
Returns a list of reciprocal lattice vectors within distance dq of a vector q.
Definition: Lattice3D.cpp:92
~Lattice3D() override
kvector_t m_b
Definition: Lattice3D.h:77
void onChange() override
Action to be taken in inherited class when a parameter has changed.
Definition: Lattice3D.cpp:50
Lattice3D()=delete
void initialize()
Initializes cached data.
Definition: Lattice3D.cpp:40
Vector transformations in three dimensions.
Definition: Transform3D.h:26
ivector_t transformed(const ivector_t &v) const
Return transformed vector v.