oILAB
Loading...
Searching...
No Matches
ReciprocalLatticeDirectionBindings.h
Go to the documentation of this file.
1//
2// Created by Nikhil Chandra Admal on 7/4/25.
3//
4
5#ifndef OILAB_RECIPROCALLATTICEDIRECTIONBINDINGS_H
6#define OILAB_RECIPROCALLATTICEDIRECTIONBINDINGS_H
7
8#include <pybind11/pybind11.h>
9#include <pybind11/operators.h>
10#include <LatticeModule.h>
11#include <pybind11/numpy.h>
12#include <pybind11/eigen.h>
13
14namespace pyoilab{
15 template<int dim>
17 {
22
23 using IntScalarType = long long int;
24 using MatrixDimD = Eigen::Matrix<double, dim, dim>;
25 using VectorDimD = Eigen::Matrix<double, dim, 1>;
26 using VectorDimI = Eigen::Matrix<IntScalarType, dim, 1>;
27 using MatrixDimI = Eigen::Matrix<IntScalarType, dim, dim>;
28 public:
30
35
40
42 return rld.cartesian();
43 }
44
48
50 return rld.dot(other.lv);
51 }
52
53 double planeSpacing() const {
54 return rld.planeSpacing();
55 }
56
57 int stacking() const {
58 return rld.stacking();
59 }
60 };
61
62 template<int dim>
63 void bind_ReciprocalLatticeDirection(py::module_ &m) {
64 using MatrixDimD = Eigen::Matrix<double, dim, dim>;
65 using VectorDimD = Eigen::Matrix<double, dim, 1>;
66 using VectorDimI = Eigen::Matrix<long long int, dim, 1>;
67
68 using Lattice = gbLAB::Lattice<dim>;
69 using ReciprocalLatticeDirection = gbLAB::ReciprocalLatticeDirection<dim>;
70 using LatticeVector = gbLAB::LatticeVector<dim>;
71
74
75 py::class_<PyReciprocalLatticeDirection>(m, ("ReciprocalLatticeDirection" + std::to_string(dim) + "D").c_str())
76 .def(py::init<const PyReciprocalLatticeVector&>())
77 .def(py::init<const PyReciprocalLatticeDirection&>())
78 .def("reciprocalLatticeVector",[](const PyReciprocalLatticeDirection& rld) {
80 },"Get the reciprocal lattice vector")
81 .def("cartesian",&PyReciprocalLatticeDirection::cartesian, "Cartesian coordinates of the reciprocal lattice direction.")
82 .def("integerCoordinates",&PyReciprocalLatticeDirection::integerCoordinates,"Integer coordinates of the reciprocal lattice direction.")
83 .def("dot",&PyReciprocalLatticeDirection::dot,"dot product with a lattice vector")
84 .def("planeSpacing",&PyReciprocalLatticeDirection::planeSpacing,"distance between two consecutive planes")
85 .def("stacking",&PyReciprocalLatticeDirection::stacking,"stacking");
86
87 }
88}
89#endif //OILAB_RECIPROCALLATTICEDIRECTIONBINDINGS_H
Lattice class.
Definition Lattice.h:34
LatticeVector class.
IntScalarType dot(const LatticeVector< dim > &other) const
Eigen::Matrix< IntScalarType, dim, dim > MatrixDimI
Eigen::Matrix< IntScalarType, dim, 1 > VectorDimI
const ReciprocalLatticeVector & reciprocalLatticeVector() const
PyReciprocalLatticeDirection(const ReciprocalLatticeDirection &rld)
PyReciprocalLatticeDirection(const PyReciprocalLatticeDirection &prld)=default
PyReciprocalLatticeDirection(const ReciprocalLatticeVector &rlv)
PyReciprocalLatticeDirection(const PyReciprocalLatticeVector &prlv)
IntScalarType dot(const PyLatticeVector< dim > &other) const
void bind_ReciprocalLatticeDirection(py::module_ &m)
double planeSpacing() const
Returns the spacing between two consecutive lattice planes.
int stacking() const
Returns the number of planes in the stacking sequence.
const ReciprocalLatticeVector< dim > & reciprocalLatticeVector() const
Returns a constant reference to the base class (ReciprocalLatticeVector)