oILAB
Loading...
Searching...
No Matches
ReciprocalLatticeVectorBindings.h
Go to the documentation of this file.
1//
2// Created by Nikhil Chandra Admal on 7/3/25.
3//
4
5#ifndef OILAB_RECIPROCALLATTICEVECTORBINDINGS_H
6#define OILAB_RECIPROCALLATTICEVECTORBINDINGS_H
7
8#include <pybind11/pybind11.h>
9#include <pybind11/operators.h>
10#include <pybind11/numpy.h>
11#include <pybind11/eigen.h>
12
13#include "../Lattices/LatticeModule.h"
14#include "PyLatticeModule.h"
15
16namespace py = pybind11;
17
18namespace pyoilab {
19// Had to wrap the ReciprocalLatticeVector<dim> class due to Pybind11 issues with classes with
20// Eigen base classes
21 template<int dim>
25
26 using IntScalarType = long long int;
27 using MatrixDimD = Eigen::Matrix<double, dim, dim>;
28 using VectorDimD = Eigen::Matrix<double, dim, 1>;
29 using VectorDimI = Eigen::Matrix<IntScalarType, dim, 1>;
30 using MatrixDimI = Eigen::Matrix<IntScalarType, dim, dim>;
31 public:
33
34 PyReciprocalLatticeVector(const Lattice &lattice) : rlv(lattice) {}
35
36 PyReciprocalLatticeVector(const VectorDimD &cartesianCoordinates, const Lattice &lattice) : rlv(
37 cartesianCoordinates,
38 lattice) {}
39
43
45
47
49 return PyReciprocalLatticeVector(rlv.operator+(other.rlv));
50 }
51
53 return PyReciprocalLatticeVector(rlv.operator-(other.rlv));
54 }
55
57 rlv.operator+=(other.rlv);
58 return *this;
59 }
60
62 rlv.operator-=(other.rlv);
63 return *this;
64 }
65
67 return PyReciprocalLatticeVector(rlv.operator*(scalar));
68 }
69
71 return rlv.cartesian();
72 }
73
75 return rlv;
76 }
77
79 return rlv.dot(other.lv);
80 }
81
85
90
94
95 template<int dm=dim>
96 typename std::enable_if<dm==2,PyLatticeDirection<dm>>::type
98 return PyLatticeDirection(rlv.cross(other.rlv));
99 }
100 template<int dm=dim>
101 typename std::enable_if<dm==3,PyLatticeDirection<dm>>::type
103 return PyLatticeDirection(rlv.cross(other.rlv));
104 }
105 };
106
107 template<int dim>
109 return L * scalar;
110 }
111
112
113 template<int dim>
114 void bind_ReciprocalLatticeVector(py::module_ &m) {
115 using Lattice = oILAB::Lattice<dim>;
116 using LatticeVector = oILAB::LatticeVector<dim>;
117 using ReciprocalLatticeVector = oILAB::ReciprocalLatticeVector<dim>;
118 using IntScalarType = long long int;
119 using MatrixDimD = Eigen::Matrix<double, dim, dim>;
120 using VectorDimD = Eigen::Matrix<double, dim, 1>;
121 using VectorDimI = Eigen::Matrix<IntScalarType, dim, 1>;
122 using MatrixDimI = Eigen::Matrix<IntScalarType, dim, dim>;
124
125 py::class_<PyReciprocalLatticeVector>(
126 m, ("ReciprocalLatticeVector" + std::to_string(dim) + "D").c_str())
127 .def(py::init<const Lattice &>())
128 .def(py::init<const VectorDimD &, const Lattice &>())
129 .def(py::init<const VectorDimI &, const Lattice &>())
130 .def(py::init<const PyReciprocalLatticeVector &>())
131 .def("cartesian", &PyReciprocalLatticeVector::cartesian)
132 .def("integerCoordinates",
134 .def(py::self + py::self)
135 .def(py::self - py::self)
136 .def(
137 "__mul__",
138 [](const PyReciprocalLatticeVector &self,
139 const IntScalarType &scalar) { return self * scalar; },
140 py::is_operator())
141 .def(
142 "__rmul__",
143 [](const PyReciprocalLatticeVector &self,
144 const IntScalarType &scalar) { return scalar * self; },
145 py::is_operator())
146 .def(py::self -= py::self)
147 .def(py::self += py::self)
149 // note that cross is a template member function
150 .def("cross", &PyReciprocalLatticeVector::template cross<dim>);
151 }
152
153}
154#endif //OILAB_RECIPROCALLATTICEVECTORBINDINGS_H
Lattice class.
Definition Lattice.h:31
LatticeVector class.
IntScalarType closestPlaneIndexOfPoint(const VectorDimD &P) const
IntScalarType dot(const LatticeVector< dim > &other) const
IntScalarType planeIndexOfPoint(const VectorDimD &P) const
std::enable_if< dm==2, LatticeDirection< dim > >::type cross(const ReciprocalLatticeVector< dim > &other) const
PyReciprocalLatticeVector operator+(const PyReciprocalLatticeVector &other) const
IntScalarType closestPlaneIndexOfPoint(const VectorDimD &p) const
PyReciprocalLatticeVector operator+=(const PyReciprocalLatticeVector &other)
IntScalarType dot(const PyLatticeVector< dim > &other) const
PyReciprocalLatticeVector operator-=(const PyReciprocalLatticeVector &other)
std::enable_if< dm==2, PyLatticeDirection< dm > >::type cross(const PyReciprocalLatticeVector< dm > &other)
PyReciprocalLatticeVector(const PyReciprocalLatticeVector &)=default
std::enable_if< dm==3, PyLatticeDirection< dm > >::type cross(const PyReciprocalLatticeVector< dm > &other)
PyReciprocalLatticeVector(const ReciprocalLatticeVector &rlv)
PyReciprocalLatticeVector operator*(const IntScalarType &scalar) const
PyReciprocalLatticeVector(const VectorDimI &integerCoordinates, const Lattice &lattice)
IntScalarType planeIndexOfPoint(const VectorDimD &p)
IntScalarType planeIndexOfPoint(const oILAB::LatticeVector< dim > &p) const
Eigen::Matrix< IntScalarType, dim, 1 > VectorDimI
PyReciprocalLatticeVector operator-(const PyReciprocalLatticeVector &other) const
PyReciprocalLatticeVector(const VectorDimD &cartesianCoordinates, const Lattice &lattice)
Eigen::Matrix< IntScalarType, dim, dim > MatrixDimI
Eigen::Matrix< double, dim, dim > MatrixDimD
void bind_ReciprocalLatticeVector(py::module_ &m)
PyLatticeVector< dim > operator*(const long long int &scalar, const PyLatticeVector< dim > &L)