BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
Lattice.h
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Sample/Lattice/Lattice.h
6 //! @brief Defines 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 
15 #ifndef BORNAGAIN_CORE_LATTICE_LATTICE_H
16 #define BORNAGAIN_CORE_LATTICE_LATTICE_H
17 
18 #include "Param/Node/INode.h"
19 #include <vector>
20 
21 class ISelectionRule;
22 class Transform3D;
23 
24 //! A lattice with three basis vectors.
25 //! @ingroup samples
26 
27 class Lattice : public INode
28 {
29 public:
30  Lattice();
31  Lattice(const kvector_t a1, const kvector_t a2, const kvector_t a3);
32  Lattice(const Lattice& lattice);
33  ~Lattice() override;
34 
35  void accept(INodeVisitor* visitor) const override { visitor->visit(this); }
36 
37  //! Creates transformed lattice
38  Lattice createTransformedLattice(const Transform3D& transform) const;
39 
40  //! Initializes cached data
41  void initialize() const;
42 
43  //! Returns basis vector a
44  kvector_t getBasisVectorA() const { return m_a; }
45 
46  //! Returns basis vector b
47  kvector_t getBasisVectorB() const { return m_b; }
48 
49  //! Returns basis vector c
50  kvector_t getBasisVectorC() const { return m_c; }
51 
52  //! Resets the basis vectors
53  void resetBasis(const kvector_t a1, const kvector_t a2, const kvector_t a3);
54 
55  //! Returns normalized direction corresponding to the given Miller indices
56  kvector_t getMillerDirection(double h, double k, double l) const;
57 
58  //! Returns the volume of the unit cell
59  double volume() const;
60 
61  //! Returns the reciprocal basis vectors
62  void getReciprocalLatticeBasis(kvector_t& b1, kvector_t& b2, kvector_t& b3) const;
63 
64  //! Returns the nearest lattice point from a given vector
66 
67  //! Returns the nearest reciprocal lattice point from a given vector
69 
70  //! Computes a list of reciprocal lattice vectors within a specified distance of a given vector
71  std::vector<kvector_t> reciprocalLatticeVectorsWithinRadius(const kvector_t input_vector,
72  double radius) const;
73 
74  //! Sets a selection rule for the reciprocal vectors
75  void setSelectionRule(const ISelectionRule& p_selection_rule);
76 
77  //! Returns a primitive cubic (cP) lattice with edge length a.
78  static Lattice createCubicLattice(double a);
79 
80  //! Returns a face-centered cubic (cF) lattice with edge length a.
81  static Lattice createFCCLattice(double a);
82 
83  //! Returns a primitive hexagonal (hP) lattice with hexagonal edge a and height c.
84  static Lattice createHexagonalLattice(double a, double c);
85 
86  //! TODO: Clarify how this is meant: HCP is not a Bravais lattice.
87  static Lattice createHCPLattice(double a, double c);
88 
89  //! Returns a primitive tetragonal (tP) lattice with square base edge a and height c.
90  static Lattice createTetragonalLattice(double a, double c);
91 
92  //! Returns a body-centered cubic (cI) lattice with edge length a.
93  //! TODO: Clarify meaning of c
94  static Lattice createBCTLattice(double a, double c);
95 
96  void onChange() override;
97 
98 private:
99  Lattice& operator=(const Lattice& lattice);
100 
101  void registerBasisVectors();
102 
103  std::vector<kvector_t> vectorsWithinRadius(const kvector_t input_vector,
104  const ivector_t& nearest_coords, double radius,
105  const kvector_t v1, const kvector_t v2,
106  const kvector_t v3, const kvector_t rec1,
107  const kvector_t rec2, const kvector_t rec3) const;
108 
109  void computeReciprocalVectors() const;
110  static void computeInverseVectors(const kvector_t v1, const kvector_t v2, const kvector_t v3,
111  kvector_t o1, kvector_t o2, kvector_t o3);
113  kvector_t m_a, m_b, m_c; //!< Basis vectors in real space
114  mutable kvector_t m_ra, m_rb, m_rc; //!< Cache of basis vectors in reciprocal space
115  //! Boolean indicating if the reciprocal vectors are already initialized in the cache
116  mutable bool m_cache_ok;
117 };
118 
119 #endif // BORNAGAIN_CORE_LATTICE_LATTICE_H
Defines class INode.
Visitor interface to visit ISample objects.
Definition: INodeVisitor.h:149
virtual void visit(const BasicLattice *)
Definition: INodeVisitor.h:154
Base class for tree-like structures containing parameterized objects.
Definition: INode.h:49
Pure virtual base class for selection rules.
A lattice with three basis vectors.
Definition: Lattice.h:28
bool m_cache_ok
Boolean indicating if the reciprocal vectors are already initialized in the cache.
Definition: Lattice.h:116
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
void resetBasis(const kvector_t a1, const kvector_t a2, const kvector_t a3)
Resets the basis vectors.
Definition: Lattice.cpp:72
void registerBasisVectors()
Definition: Lattice.cpp:190
ivector_t getNearestLatticeVectorCoordinates(const kvector_t vector_in) const
Returns the nearest lattice point from a given vector.
Definition: Lattice.cpp:104
static void computeInverseVectors(const kvector_t v1, const kvector_t v2, const kvector_t v3, kvector_t o1, kvector_t o2, kvector_t o3)
Definition: Lattice.cpp:237
ivector_t getNearestReciprocalLatticeVectorCoordinates(const kvector_t vector_in) const
Returns the nearest reciprocal lattice point from a given vector.
Definition: Lattice.cpp:115
kvector_t m_b
Definition: Lattice.h:113
void setSelectionRule(const ISelectionRule &p_selection_rule)
Sets a selection rule for the reciprocal vectors.
Definition: Lattice.cpp:277
Lattice()
Definition: Lattice.cpp:22
kvector_t getMillerDirection(double h, double k, double l) const
Returns normalized direction corresponding to the given Miller indices.
Definition: Lattice.cpp:80
static Lattice createHCPLattice(double a, double c)
TODO: Clarify how this is meant: HCP is not a Bravais lattice.
Definition: Lattice.cpp:161
kvector_t m_c
Basis vectors in real space.
Definition: Lattice.h:113
void initialize() const
Initializes cached data.
Definition: Lattice.cpp:66
kvector_t getBasisVectorB() const
Returns basis vector b.
Definition: Lattice.h:47
static Lattice createBCTLattice(double a, double c)
Returns a body-centered cubic (cI) lattice with edge length a.
Definition: Lattice.cpp:177
std::vector< kvector_t > 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: Lattice.cpp:126
kvector_t m_rc
Cache of basis vectors in reciprocal space.
Definition: Lattice.h:114
kvector_t getBasisVectorC() const
Returns basis vector c.
Definition: Lattice.h:50
void accept(INodeVisitor *visitor) const override
Calls the INodeVisitor's visit method.
Definition: Lattice.h:35
Lattice createTransformedLattice(const Transform3D &transform) const
Creates transformed lattice.
Definition: Lattice.cpp:55
Lattice & operator=(const Lattice &lattice)
ISelectionRule * mp_selection_rule
Definition: Lattice.h:112
double volume() const
Returns the volume of the unit cell.
Definition: Lattice.cpp:88
static Lattice createTetragonalLattice(double a, double c)
Returns a primitive tetragonal (tP) lattice with square base edge a and height c.
Definition: Lattice.cpp:169
void computeReciprocalVectors() const
Definition: Lattice.cpp:199
kvector_t getBasisVectorA() const
Returns basis vector a.
Definition: Lattice.h:44
static Lattice createFCCLattice(double a)
Returns a face-centered cubic (cF) lattice with edge length a.
Definition: Lattice.cpp:144
kvector_t m_ra
Definition: Lattice.h:114
static Lattice createHexagonalLattice(double a, double c)
Returns a primitive hexagonal (hP) lattice with hexagonal edge a and height c.
Definition: Lattice.cpp:153
static Lattice createCubicLattice(double a)
Returns a primitive cubic (cP) lattice with edge length a.
Definition: Lattice.cpp:136
kvector_t m_a
Definition: Lattice.h:113
~Lattice() override
Definition: Lattice.cpp:50
void getReciprocalLatticeBasis(kvector_t &b1, kvector_t &b2, kvector_t &b3) const
Returns the reciprocal basis vectors.
Definition: Lattice.cpp:93
void onChange() override
Action to be taken in inherited class when a parameter has changed.
Definition: Lattice.cpp:185
kvector_t m_rb
Definition: Lattice.h:114
Vector transformations in three dimensions.
Definition: Transform3D.h:28
const double radius(5 *Units::nanometer)