BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
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"
17 #include "Base/Vector/RotMatrix.h"
19 #include <gsl/gsl_linalg.h>
20 
21 Lattice3D::Lattice3D(const R3 a, const R3 b, const R3 c)
22  : m_a(a)
23  , m_b(b)
24  , m_c(c)
25 {
27 }
28 
30  : Lattice3D(lattice.m_a, lattice.m_b, lattice.m_c)
31 {
32  if (lattice.m_selection_rule)
34 }
35 
36 Lattice3D::~Lattice3D() = default;
37 
38 Lattice3D Lattice3D::rotated(const RotMatrix& rotMatrix) const
39 {
40  R3 q1 = rotMatrix.transformed(m_a);
41  R3 q2 = rotMatrix.transformed(m_b);
42  R3 q3 = rotMatrix.transformed(m_c);
43  Lattice3D result = {q1, q2, q3};
44  if (m_selection_rule)
46  return result;
47 }
48 
49 //! Currently unused but may be useful for checks
50 R3 Lattice3D::getMillerDirection(double h, double k, double l) const
51 {
52  R3 direction = h * m_ra + k * m_rb + l * m_rc;
53  return direction.unit();
54 }
55 
57 {
58  return std::abs(m_a.dot(m_b.cross(m_c)));
59 }
60 
61 //! Currently only used in tests
62 void Lattice3D::reciprocalLatticeBasis(R3& ra, R3& rb, R3& rc) const
63 {
64  ra = m_ra;
65  rb = m_rb;
66  rc = m_rc;
67 }
68 
70 {
71  return {(int)std::lround(q.dot(m_a) / M_TWOPI), (int)std::lround(q.dot(m_b) / M_TWOPI),
72  (int)std::lround(q.dot(m_c) / M_TWOPI)};
73 }
74 
75 std::vector<R3> Lattice3D::reciprocalLatticeVectorsWithinRadius(const R3 q, double dq) const
76 {
77  I3 nearest_coords = nearestReciprocalLatticeVectorCoordinates(q);
78 
79  int max_X = std::lround(m_a.mag() * dq / M_TWOPI);
80  int max_Y = std::lround(m_b.mag() * dq / M_TWOPI);
81  int max_Z = std::lround(m_c.mag() * dq / M_TWOPI);
82 
83  std::vector<R3> result;
84  for (int index_X = -max_X; index_X <= max_X; ++index_X) {
85  for (int index_Y = -max_Y; index_Y <= max_Y; ++index_Y) {
86  for (int index_Z = -max_Z; index_Z <= max_Z; ++index_Z) {
87  I3 coords = I3(index_X, index_Y, index_Z) + nearest_coords;
88  if (m_selection_rule && !m_selection_rule->coordinateSelected(coords))
89  continue;
90  R3 latticePoint = coords.x() * m_ra + coords.y() * m_rb + coords.z() * m_rc;
91  if ((latticePoint - q).mag() <= dq)
92  result.push_back(latticePoint);
93  }
94  }
95  }
96  return result;
97 }
98 
100 {
101  R3 q23 = m_b.cross(m_c);
102  R3 q31 = m_c.cross(m_a);
103  R3 q12 = m_a.cross(m_b);
104  m_ra = M_TWOPI / m_a.dot(q23) * q23;
105  m_rb = M_TWOPI / m_b.dot(q31) * q31;
106  m_rc = M_TWOPI / m_c.dot(q12) * q12;
107 }
108 
109 void Lattice3D::setSelectionRule(const ISelectionRule& selection_rule)
110 {
111  m_selection_rule.reset(selection_rule.clone());
112 }
Defines M_PI and some more mathematical constants.
#define M_TWOPI
Definition: Constants.h:54
Defines classes ISelectionRule, SimpleSelectionRule.
Defines class Lattice.
Declares class RotMatrix.
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:30
I3 nearestReciprocalLatticeVectorCoordinates(R3 q) const
Returns the nearest reciprocal lattice point from a given vector.
Definition: Lattice3D.cpp:69
Lattice3D rotated(const RotMatrix &rotMatrix) const
Creates rotated lattice.
Definition: Lattice3D.cpp:38
double unitCellVolume() const
Returns the volume of the unit cell.
Definition: Lattice3D.cpp:56
R3 m_rc
Cache of basis vectors in reciprocal space.
Definition: Lattice3D.h:75
R3 getMillerDirection(double h, double k, double l) const
Returns normalized direction corresponding to the given Miller indices.
Definition: Lattice3D.cpp:50
void computeReciprocalVectors() const
Definition: Lattice3D.cpp:99
void setSelectionRule(const ISelectionRule &selection_rule)
Sets a selection rule for the reciprocal vectors.
Definition: Lattice3D.cpp:109
R3 m_c
Basis vectors in real space.
Definition: Lattice3D.h:72
std::unique_ptr< ISelectionRule > m_selection_rule
Definition: Lattice3D.h:73
~Lattice3D() override
void reciprocalLatticeBasis(R3 &ra, R3 &rb, R3 &rc) const
Returns the reciprocal basis vectors.
Definition: Lattice3D.cpp:62
std::vector< R3 > reciprocalLatticeVectorsWithinRadius(R3 q, double dq) const
Returns a list of reciprocal lattice vectors within distance dq of a vector q.
Definition: Lattice3D.cpp:75
Lattice3D(R3 a, R3 b, R3 c)
Definition: Lattice3D.cpp:21
Rotation matrix in three dimensions. Represents group SO(3). Internal parameterization based on quate...
Definition: RotMatrix.h:25
T transformed(const T &v) const
Return transformed vector v.
Definition: RotMatrix.cpp:76