BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
Lattice3D.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Lattice/Lattice3D.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_SAMPLE_LATTICE_LATTICE3D_H
16 #define BORNAGAIN_SAMPLE_LATTICE_LATTICE3D_H
17 
18 #include "Param/Node/INode.h"
19 #include <memory>
20 #include <vector>
21 
22 class ISelectionRule;
23 class Transform3D;
24 
25 //! A Bravais lattice, characterized by three basis vectors, and optionally an ISelectionRule.
26 
27 //! @ingroup samples
28 
29 class Lattice3D : public INode {
30 public:
31  Lattice3D() = delete;
32  Lattice3D(const kvector_t a, const kvector_t b, const kvector_t c);
33  Lattice3D(const Lattice3D& lattice);
34  ~Lattice3D() override;
35  Lattice3D& operator=(const Lattice3D&) = delete;
36 
37  void accept(INodeVisitor* visitor) const override { visitor->visit(this); }
38 
39  //! Creates transformed lattice
40  Lattice3D transformed(const Transform3D& transform) const;
41 
42  //! Initializes cached data
43  void initialize();
44 
45  //! Returns basis vector a
46  kvector_t getBasisVectorA() const { return m_a; }
47 
48  //! Returns basis vector b
49  kvector_t getBasisVectorB() const { return m_b; }
50 
51  //! Returns basis vector c
52  kvector_t getBasisVectorC() const { return m_c; }
53 
54  //! Returns normalized direction corresponding to the given Miller indices
55  kvector_t getMillerDirection(double h, double k, double l) const;
56 
57  //! Returns the volume of the unit cell
58  double unitCellVolume() const;
59 
60  //! Returns the reciprocal basis vectors
61  void getReciprocalLatticeBasis(kvector_t& ra, kvector_t& rb, kvector_t& rc) const;
62 
63  //! Returns the nearest reciprocal lattice point from a given vector
65 
66  //! Returns a list of reciprocal lattice vectors within distance dq of a vector q
67  std::vector<kvector_t> reciprocalLatticeVectorsWithinRadius(const kvector_t q, double dq) const;
68 
69  //! Sets a selection rule for the reciprocal vectors
70  void setSelectionRule(const ISelectionRule& selection_rule);
71 
72 private:
73  void onChange() override;
74 
75  void computeReciprocalVectors() const;
76 
77  kvector_t m_a, m_b, m_c; //!< Basis vectors in real space
78  std::unique_ptr<ISelectionRule> m_selection_rule;
79 
80  mutable kvector_t m_ra, m_rb, m_rc; //!< Cache of basis vectors in reciprocal space
81 };
82 
83 #endif // BORNAGAIN_SAMPLE_LATTICE_LATTICE3D_H
Defines interface INode.
Visitor interface to visit ISampleNode objects.
Definition: INodeVisitor.h:146
virtual void visit(const BasicLattice2D *)
Definition: INodeVisitor.h:151
Base class for tree-like structures containing parameterized objects.
Definition: INode.h:49
Abstract base class for selection rules.
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 getBasisVectorC() const
Returns basis vector c.
Definition: Lattice3D.h:52
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
kvector_t getBasisVectorA() const
Returns basis vector a.
Definition: Lattice3D.h:46
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 & operator=(const Lattice3D &)=delete
~Lattice3D() override
kvector_t getBasisVectorB() const
Returns basis vector b.
Definition: Lattice3D.h:49
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
void accept(INodeVisitor *visitor) const override
Calls the INodeVisitor's visit method.
Definition: Lattice3D.h:37
Lattice3D()=delete
void initialize()
Initializes cached data.
Definition: Lattice3D.cpp:40
Vector transformations in three dimensions.
Definition: Transform3D.h:26