Line data Source code
1 : /*
2 : * Copyright (c) 2020 Sebastian Weber, Henri Menke, Alexander Papageorge. All rights reserved.
3 : *
4 : * This file is part of the pairinteraction library.
5 : *
6 : * The pairinteraction library is free software: you can redistribute it and/or modify
7 : * it under the terms of the GNU Lesser General Public License as published by
8 : * the Free Software Foundation, either version 3 of the License, or
9 : * (at your option) any later version.
10 : *
11 : * The pairinteraction library is distributed in the hope that it will be useful,
12 : * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 : * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 : * GNU Lesser General Public License for more details.
15 : *
16 : * You should have received a copy of the GNU Lesser General Public License
17 : * along with the pairinteraction library. If not, see <http://www.gnu.org/licenses/>.
18 : */
19 :
20 : #include "EigenCompat.hpp"
21 : #include <Eigen/SparseCore>
22 :
23 : #include "Constants.hpp"
24 : #include "MatrixElementCache.hpp"
25 : #include "State.hpp"
26 : #include "SystemBase.hpp"
27 : #include "SystemOne.hpp"
28 : #include "SystemTwo.hpp"
29 : #include <jlcxx/const_array.hpp>
30 : #include <jlcxx/jlcxx.hpp>
31 :
32 : namespace jlcxx {
33 : template <>
34 : struct IsBits<method_t> : std::true_type {};
35 : template <>
36 : struct IsBits<parity_t> : std::true_type {};
37 :
38 4 : jlcxx::Array<double> get_array_from_evd_t(Eigen::VectorX<double> overlap) {
39 4 : jlcxx::Array<double> ret;
40 788 : for (unsigned i = 0; i < overlap.size(); i++) {
41 784 : ret.push_back(overlap.data()[i]);
42 : }
43 4 : return ret;
44 : }
45 :
46 : template <typename T>
47 : struct BuildParameterList<Eigen::SparseMatrix<T>> {
48 : typedef ParameterList<T> type;
49 : };
50 : } // namespace jlcxx
51 :
52 6 : JLCXX_MODULE define_julia_module(jlcxx::Module &pi) {
53 6 : pi.add_bits<method_t>("method_t");
54 6 : pi.set_const("NUMEROV", NUMEROV);
55 6 : pi.set_const("WHITTAKER", WHITTAKER);
56 :
57 6 : pi.add_bits<parity_t>("parity_t");
58 6 : pi.set_const("NA", NA);
59 6 : pi.set_const("EVEN", EVEN);
60 6 : pi.set_const("ODD", ODD);
61 :
62 6 : pi.set_const("ARB", ARB);
63 :
64 12 : pi.add_type<MatrixElementCache>("MatrixElementCache")
65 6 : .constructor<std::string>()
66 12 : .method("getElectricDipole", &MatrixElementCache::getElectricDipole)
67 : .method(
68 : "getElectricMultipole",
69 : static_cast<double (MatrixElementCache::*)(const StateOne &, const StateOne &, int)>(
70 12 : &MatrixElementCache::getElectricMultipole))
71 18 : .method("getDiamagnetism", &MatrixElementCache::getDiamagnetism)
72 12 : .method("getMagneticDipole", &MatrixElementCache::getMagneticDipole)
73 : .method("getElectricMultipole",
74 : static_cast<double (MatrixElementCache::*)(const StateOne &, const StateOne &, int,
75 : int)>(
76 12 : &MatrixElementCache::getElectricMultipole))
77 18 : .method("getRadial", &MatrixElementCache::getRadial)
78 : .method("precalculateElectricMomentum",
79 0 : [](MatrixElementCache &mec, jlcxx::ArrayRef<jl_value_t *> basis_one_jl, int q) {
80 0 : std::vector<StateOne> basis_one;
81 0 : for (unsigned i = 0; i < basis_one_jl.size(); i++) {
82 0 : const StateOne s = *jlcxx::unbox_wrapped_ptr<StateOne>(basis_one_jl[i]);
83 0 : basis_one.push_back(s);
84 : }
85 0 : const std::vector<StateOne> &basis_one_const = basis_one;
86 0 : mec.precalculateElectricMomentum(basis_one_const, q);
87 12 : })
88 : .method("precalculateMagneticMomentum",
89 0 : [](MatrixElementCache &mec, jlcxx::ArrayRef<jl_value_t *> basis_one_jl, int q) {
90 0 : std::vector<StateOne> basis_one;
91 0 : for (unsigned i = 0; i < basis_one_jl.size(); i++) {
92 0 : const StateOne s = *jlcxx::unbox_wrapped_ptr<StateOne>(basis_one_jl[i]);
93 0 : basis_one.push_back(s);
94 : }
95 0 : const std::vector<StateOne> &basis_one_const = basis_one;
96 0 : mec.precalculateMagneticMomentum(basis_one_const, q);
97 12 : })
98 : .method(
99 : "precalculateDiamagnetism",
100 0 : [](MatrixElementCache &mec, jlcxx::ArrayRef<jl_value_t *> basis_one_jl, int k, int q) {
101 0 : std::vector<StateOne> basis_one;
102 0 : for (unsigned i = 0; i < basis_one_jl.size(); i++) {
103 0 : const StateOne s = *jlcxx::unbox_wrapped_ptr<StateOne>(basis_one_jl[i]);
104 0 : basis_one.push_back(s);
105 : }
106 0 : const std::vector<StateOne> &basis_one_const = basis_one;
107 0 : mec.precalculateDiamagnetism(basis_one_const, k, q);
108 12 : })
109 : .method("precalculateMultipole",
110 0 : [](MatrixElementCache &mec, jlcxx::ArrayRef<jl_value_t *> basis_one_jl, int k) {
111 0 : std::vector<StateOne> basis_one;
112 0 : for (unsigned i = 0; i < basis_one_jl.size(); i++) {
113 0 : const StateOne s = *jlcxx::unbox_wrapped_ptr<StateOne>(basis_one_jl[i]);
114 0 : basis_one.push_back(s);
115 : }
116 0 : const std::vector<StateOne> &basis_one_const = basis_one;
117 0 : mec.precalculateMultipole(basis_one_const, k);
118 12 : })
119 : .method("precalculateRadial",
120 0 : [](MatrixElementCache &mec, jlcxx::ArrayRef<jl_value_t *> basis_one_jl, int k) {
121 0 : std::vector<StateOne> basis_one;
122 0 : for (unsigned i = 0; i < basis_one_jl.size(); i++) {
123 0 : const StateOne s = *jlcxx::unbox_wrapped_ptr<StateOne>(basis_one_jl[i]);
124 0 : basis_one.push_back(s);
125 : }
126 0 : const std::vector<StateOne> &basis_one_const = basis_one;
127 0 : mec.precalculateRadial(basis_one_const, k);
128 12 : })
129 12 : .method("setDefectDB", &MatrixElementCache::setDefectDB)
130 12 : .method("setMethod", &MatrixElementCache::setMethod)
131 12 : .method("loadElectricDipoleDB", &MatrixElementCache::loadElectricDipoleDB)
132 6 : .method("size", &MatrixElementCache::size);
133 :
134 12 : pi.add_type<StateOne>("StateOne")
135 6 : .constructor<std::string, int, int, float, float>()
136 6 : .constructor<std::string>()
137 12 : .method("string", static_cast<std::string (StateOne::*)() const>(&StateOne::str))
138 18 : .method("getN", &StateOne::getN)
139 12 : .method("getL", &StateOne::getL)
140 12 : .method("getJ", &StateOne::getJ)
141 12 : .method("getM", &StateOne::getM)
142 12 : .method("getS", &StateOne::getS)
143 12 : .method("getSpecies", &StateOne::getSpecies)
144 12 : .method("getElement", &StateOne::getElement)
145 12 : .method("getEnergy", static_cast<double (StateOne::*)() const>(&StateOne::getEnergy))
146 18 : .method("getNStar", static_cast<double (StateOne::*)() const>(&StateOne::getNStar))
147 18 : .method("getLabel", &StateOne::getLabel)
148 12 : .method("isArtificial", &StateOne::isArtificial)
149 12 : .method("isGeneralized", &StateOne::isGeneralized)
150 12 : .method("getHash", &StateOne::getHash)
151 12 : .method("getReflected", &StateOne::getReflected)
152 12 : .method("==", &StateOne::operator==)
153 12 : .method("^", &StateOne::operator^)
154 12 : .method("!=", &StateOne::operator!=)
155 12 : .method("<", &StateOne::operator<)
156 6 : .method("<=", &StateOne::operator<=);
157 :
158 12 : pi.add_type<StateTwo>("StateTwo")
159 6 : .constructor<StateOne, StateOne>()
160 : .method("StateTwo",
161 1 : [](jlcxx::ArrayRef<std::string> init_arr) {
162 2 : std::array<std::string, 2> str_arr = {init_arr[0], init_arr[1]};
163 1 : StateTwo s = StateTwo(str_arr);
164 2 : return s;
165 12 : })
166 : .method("StateTwo",
167 0 : [](jlcxx::ArrayRef<std::string> species, jlcxx::ArrayRef<int> ns,
168 : jlcxx::ArrayRef<int> ls, jlcxx::ArrayRef<float> js, jlcxx::ArrayRef<float> ms) {
169 0 : std::array<std::string, 2> species_arr = {species[0], species[1]};
170 0 : std::array<int, 2> n_arr = {ns[0], ns[1]};
171 0 : std::array<int, 2> l_arr = {ls[0], ls[1]};
172 0 : std::array<float, 2> j_arr = {js[0], js[1]};
173 0 : std::array<float, 2> m_arr = {ms[0], ms[1]};
174 0 : StateTwo s = StateTwo(species_arr, n_arr, l_arr, j_arr, m_arr);
175 0 : return s;
176 12 : })
177 12 : .method("str", &StateTwo::str)
178 : .method("getN",
179 1 : [](StateTwo &s) {
180 1 : std::array<int, 2> n_arr = s.getN();
181 2 : return std::make_tuple(n_arr[0], n_arr[1]);
182 12 : })
183 12 : .method("getN", static_cast<const int &(StateTwo::*)(int) const>(&StateTwo::getN))
184 : .method("getL",
185 1 : [](StateTwo &s) {
186 1 : std::array<int, 2> l_arr = s.getL();
187 2 : return std::make_tuple(l_arr[0], l_arr[1]);
188 18 : })
189 12 : .method("getL", static_cast<const int &(StateTwo::*)(int) const>(&StateTwo::getL))
190 : .method("getJ",
191 1 : [](StateTwo &s) {
192 1 : std::array<float, 2> j_arr = s.getJ();
193 2 : return std::make_tuple(j_arr[0], j_arr[1]);
194 18 : })
195 12 : .method("getJ", static_cast<const float &(StateTwo::*)(int) const>(&StateTwo::getJ))
196 : .method("getM",
197 1 : [](StateTwo &s) {
198 1 : std::array<float, 2> m_arr = s.getM();
199 2 : return std::make_tuple(m_arr[0], m_arr[1]);
200 18 : })
201 12 : .method("getM", static_cast<const float &(StateTwo::*)(int) const>(&StateTwo::getM))
202 : .method("getS",
203 1 : [](StateTwo &s) {
204 1 : std::array<float, 2> s_arr = s.getS();
205 2 : return std::make_tuple(s_arr[0], s_arr[1]);
206 18 : })
207 12 : .method("getS", static_cast<const float &(StateTwo::*)(int) const>(&StateTwo::getS))
208 : .method("getSpecies",
209 1 : [](StateTwo &s) {
210 2 : std::array<std::string, 2> species_arr = s.getSpecies();
211 2 : return std::make_tuple(species_arr[0], species_arr[1]);
212 18 : })
213 : .method("getSpecies",
214 12 : static_cast<const std::string &(StateTwo::*)(int) const>(&StateTwo::getSpecies))
215 : .method("getElement",
216 1 : [](StateTwo &s) {
217 2 : std::array<std::string, 2> el_arr = s.getElement();
218 2 : return std::make_tuple(el_arr[0], el_arr[1]);
219 18 : })
220 : .method("getElement",
221 12 : static_cast<const std::string &(StateTwo::*)(int) const>(&StateTwo::getElement))
222 18 : .method("getEnergy", static_cast<double (StateTwo::*)() const>(&StateTwo::getEnergy))
223 18 : .method("getEnergy", static_cast<double (StateTwo::*)(int) const>(&StateTwo::getEnergy))
224 : .method("getNStar",
225 0 : [](StateTwo &s) {
226 0 : std::array<double, 2> ns_arr = s.getNStar();
227 0 : return std::make_tuple(ns_arr[0], ns_arr[1]);
228 18 : })
229 12 : .method("getNStar", static_cast<double (StateTwo::*)(int) const>(&StateTwo::getNStar))
230 18 : .method("getLeRoyRadius", &StateTwo::getLeRoyRadius)
231 : .method("getLabel",
232 1 : [](StateTwo &s) {
233 2 : std::array<std::string, 2> l_arr = s.getLabel();
234 2 : return std::make_tuple(l_arr[0], l_arr[1]);
235 12 : })
236 : .method("getLabel",
237 12 : static_cast<const std::string &(StateTwo::*)(int) const>(&StateTwo::getLabel))
238 : .method("isArtificial",
239 1 : [](StateTwo &s) {
240 1 : std::array<bool, 2> b_arr = s.isArtificial();
241 2 : return std::make_tuple(b_arr[0], b_arr[1]);
242 18 : })
243 12 : .method("isArtificial", static_cast<bool (StateTwo::*)(int) const>(&StateTwo::isArtificial))
244 : .method("isGeneralized",
245 0 : [](StateTwo &s) {
246 0 : std::array<bool, 2> b_arr = s.isGeneralized();
247 0 : return std::make_tuple(b_arr[0], b_arr[1]);
248 18 : })
249 : .method("isGeneralized",
250 12 : static_cast<bool (StateTwo::*)(int) const>(&StateTwo::isGeneralized))
251 18 : .method("getFirstState", &StateTwo::getFirstState)
252 12 : .method("getSecondState", &StateTwo::getSecondState)
253 12 : .method("getHash", &StateTwo::getHash)
254 12 : .method("getReflected", &StateTwo::getReflected)
255 12 : .method("==", &StateTwo::operator==)
256 12 : .method("^", &StateTwo::operator^)
257 12 : .method("!=", &StateTwo::operator!=)
258 12 : .method("<", &StateTwo::operator<)
259 6 : .method("<=", &StateTwo::operator<=);
260 :
261 12 : pi.add_type<jlcxx::Parametric<jlcxx::TypeVar<1>, jlcxx::TypeVar<2>>>("SystemBase")
262 : .apply<SystemBase<double, StateOne>, SystemBase<std::complex<double>, StateOne>>(
263 18 : [](auto wrapped) { typedef typename decltype(wrapped)::type WrappedT; })
264 : .apply<SystemBase<double, StateTwo>, SystemBase<std::complex<double>, StateTwo>>(
265 18 : [](auto wrapped) { typedef typename decltype(wrapped)::type WrappedT; });
266 :
267 12 : pi.add_type<jlcxx::Parametric<jlcxx::TypeVar<1>>>("eigen_sparse_t")
268 : .apply<Eigen::SparseMatrix<double>, Eigen::SparseMatrix<std::complex<double>>>(
269 12 : [](auto wrapped) {
270 : typedef typename decltype(wrapped)::type WrappedT;
271 : typedef typename WrappedT::Scalar Scalar;
272 12 : constexpr bool const is_complex_v = utils::is_complex<Scalar>::value;
273 :
274 19 : wrapped.method("nonzerorealvalues", [](WrappedT &e) {
275 7 : jlcxx::Array<double> ret;
276 : if constexpr (is_complex_v) {
277 0 : for (int i = 0; i < e.nonZeros(); i++) {
278 0 : ret.push_back((e.valuePtr()[i]).real());
279 : }
280 : } else {
281 1683 : for (int i = 0; i < e.nonZeros(); i++) {
282 1676 : ret.push_back(e.valuePtr()[i]);
283 : }
284 : }
285 7 : return ret;
286 : });
287 19 : wrapped.method("nonzeroimagvalues", [](WrappedT &e) {
288 7 : jlcxx::Array<double> ret;
289 : if constexpr (is_complex_v) {
290 0 : for (int i = 0; i < e.nonZeros(); i++) {
291 0 : ret.push_back((e.valuePtr()[i]).imag());
292 : }
293 : } else {
294 1683 : for (int i = 0; i < e.nonZeros(); i++) {
295 1676 : ret.push_back(0);
296 : }
297 : }
298 7 : return ret;
299 : });
300 19 : wrapped.method("outerIndex", [](WrappedT &e) {
301 7 : jlcxx::Array<int> ret;
302 939 : for (int i = 0; i < e.outerSize(); i++) {
303 932 : ret.push_back(e.outerIndexPtr()[i]);
304 : }
305 7 : return ret;
306 : });
307 19 : wrapped.method("innerIndex", [](WrappedT &e) {
308 7 : jlcxx::Array<int> ret;
309 1683 : for (int i = 0; i < e.nonZeros(); i++) {
310 1676 : ret.push_back(e.innerIndexPtr()[i]);
311 : }
312 7 : return ret;
313 : });
314 18 : });
315 :
316 12 : pi.add_type<jlcxx::Parametric<jlcxx::TypeVar<1>>>("SystemOne")
317 12 : .apply<SystemOne<double>, SystemOne<std::complex<double>>>([](auto wrapped) {
318 : typedef typename decltype(wrapped)::type WrappedT;
319 : typedef typename WrappedT::Scalar Scalar;
320 :
321 12 : wrapped.template constructor<std::string, MatrixElementCache &>();
322 12 : wrapped.template constructor<std::string, MatrixElementCache &, bool>();
323 : // SystemBase methods
324 : // ///////////////////////////////////////////////////////////////////////
325 12 : wrapped.method("restrictEnergy", &WrappedT::restrictEnergy);
326 12 : wrapped.method("restrictN",
327 : static_cast<void (WrappedT::*)(int, int)>(&WrappedT::restrictN));
328 12 : wrapped.method("restrictN", [](WrappedT &s, jlcxx::ArrayRef<int> n_jl) {
329 0 : std::set<int> n = {n_jl[0], n_jl[1]};
330 0 : s.restrictN(n);
331 : });
332 12 : wrapped.method("restrictL",
333 : static_cast<void (WrappedT::*)(int, int)>(&WrappedT::restrictL));
334 12 : wrapped.method("restrictL", [](WrappedT &s, jlcxx::ArrayRef<int> l_jl) {
335 0 : std::set<int> l = {l_jl[0], l_jl[1]};
336 0 : s.restrictL(l);
337 : });
338 12 : wrapped.method("restrictJ",
339 : static_cast<void (WrappedT::*)(float, float)>(&WrappedT::restrictJ));
340 12 : wrapped.method("restrictJ", [](WrappedT &s, jlcxx::ArrayRef<float> j_jl) {
341 0 : std::set<float> j = {j_jl[0], j_jl[1]};
342 0 : s.restrictJ(j);
343 : });
344 12 : wrapped.method("restrictM",
345 : static_cast<void (WrappedT::*)(float, float)>(&WrappedT::restrictM));
346 12 : wrapped.method("restrictM", [](WrappedT &s, jlcxx::ArrayRef<float> m_jl) {
347 0 : std::set<float> m = {m_jl[0], m_jl[1]};
348 0 : s.restrictM(m);
349 : });
350 12 : wrapped.method("diagonalize",
351 : static_cast<void (WrappedT::*)()>(&WrappedT::diagonalize));
352 12 : wrapped.method("getHamiltonian",
353 : static_cast<Eigen::SparseMatrix<Scalar> &(WrappedT::*)()>(
354 : &WrappedT::getHamiltonian));
355 : /////////////////////////////////////////////////////////////////////////////////////////////
356 12 : wrapped.method("getSpecies", &WrappedT::getSpecies);
357 13 : wrapped.method("setEfield", [](WrappedT &s, jlcxx::ArrayRef<double> efield) {
358 1 : std::array<double, 3> field = {efield[0], efield[1], efield[2]};
359 1 : s.setEfield(field);
360 : });
361 12 : wrapped.method("setEfield",
362 0 : [](WrappedT &s, jlcxx::ArrayRef<double> efield,
363 : jlcxx::ArrayRef<double> z_axis, jlcxx::ArrayRef<double> y_axis) {
364 0 : std::array<double, 3> field = {efield[0], efield[1], efield[2]};
365 0 : std::array<double, 3> to_z_axis = {z_axis[0], z_axis[1], z_axis[2]};
366 0 : std::array<double, 3> to_y_axis = {y_axis[0], y_axis[1], y_axis[2]};
367 0 : s.setEfield(field, to_z_axis, to_y_axis);
368 : });
369 12 : wrapped.method("setEfield",
370 0 : [](WrappedT &s, jlcxx::ArrayRef<double> efield, double alpha,
371 : double beta, double gamma) {
372 0 : std::array<double, 3> field = {efield[0], efield[1], efield[2]};
373 0 : s.setEfield(field, alpha, beta, gamma);
374 : });
375 13 : wrapped.method("setBfield", [](WrappedT &s, jlcxx::ArrayRef<double> bfield) {
376 1 : std::array<double, 3> field = {bfield[0], bfield[1], bfield[2]};
377 1 : s.setBfield(field);
378 : });
379 12 : wrapped.method("setBfield",
380 0 : [](WrappedT &s, jlcxx::ArrayRef<double> bfield,
381 : jlcxx::ArrayRef<double> z_axis, jlcxx::ArrayRef<double> y_axis) {
382 0 : std::array<double, 3> field = {bfield[0], bfield[1], bfield[2]};
383 0 : std::array<double, 3> to_z_axis = {z_axis[0], z_axis[1], z_axis[2]};
384 0 : std::array<double, 3> to_y_axis = {y_axis[0], y_axis[1], y_axis[2]};
385 0 : s.setBfield(field, to_z_axis, to_y_axis);
386 : });
387 12 : wrapped.method("setBfield",
388 0 : [](WrappedT &s, jlcxx::ArrayRef<double> bfield, double alpha,
389 : double beta, double gamma) {
390 0 : std::array<double, 3> field = {bfield[0], bfield[1], bfield[2]};
391 0 : s.setBfield(field, alpha, beta, gamma);
392 : });
393 12 : wrapped.method("enableDiamagnetism", &WrappedT::enableDiamagnetism);
394 12 : wrapped.method("setConservedParityUnderReflection",
395 : &WrappedT::setConservedParityUnderReflection);
396 12 : wrapped.method("setConservedMomentaUnderRotation",
397 1 : [](WrappedT &s, jlcxx::ArrayRef<float> momenta_jl) {
398 2 : std::set<float> momenta;
399 2 : for (unsigned i = 0; i < momenta_jl.size(); i++) {
400 1 : momenta.insert(momenta_jl[i]);
401 : }
402 1 : const std::set<float> &momenta_const = momenta;
403 1 : s.setConservedMomentaUnderRotation(momenta_const);
404 : });
405 18 : });
406 :
407 12 : pi.add_type<jlcxx::Parametric<jlcxx::TypeVar<1>>>("SystemTwo")
408 12 : .apply<SystemTwo<double>, SystemTwo<std::complex<double>>>([](auto wrapped) {
409 : typedef typename decltype(wrapped)::type WrappedT;
410 : typedef typename WrappedT::Scalar Scalar;
411 :
412 : wrapped
413 12 : .template constructor<SystemOne<Scalar>, SystemOne<Scalar>, MatrixElementCache &>();
414 : wrapped.template constructor<SystemOne<Scalar>, SystemOne<Scalar>, MatrixElementCache &,
415 12 : bool>();
416 12 : wrapped.template constructor<WrappedT>();
417 : // SystemBase methods
418 : // ///////////////////////////////////////////////////////////////////////
419 12 : wrapped.method("restrictEnergy", &WrappedT::restrictEnergy);
420 12 : wrapped.method("restrictN",
421 : static_cast<void (WrappedT::*)(int, int)>(&WrappedT::restrictN));
422 12 : wrapped.method("restrictN", [](WrappedT &s, jlcxx::ArrayRef<int> n_jl) {
423 0 : std::set<int> n = {n_jl[0], n_jl[1]};
424 0 : s.restrictN(n);
425 : });
426 12 : wrapped.method("restrictL",
427 : static_cast<void (WrappedT::*)(int, int)>(&WrappedT::restrictL));
428 12 : wrapped.method("restrictL", [](WrappedT &s, jlcxx::ArrayRef<int> l_jl) {
429 0 : std::set<int> l = {l_jl[0], l_jl[1]};
430 0 : s.restrictL(l);
431 : });
432 12 : wrapped.method("restrictJ",
433 : static_cast<void (WrappedT::*)(float, float)>(&WrappedT::restrictJ));
434 12 : wrapped.method("restrictJ", [](WrappedT &s, jlcxx::ArrayRef<float> j_jl) {
435 0 : std::set<float> j = {j_jl[0], j_jl[1]};
436 0 : s.restrictJ(j);
437 : });
438 12 : wrapped.method("restrictM",
439 : static_cast<void (WrappedT::*)(float, float)>(&WrappedT::restrictM));
440 12 : wrapped.method("restrictM", [](WrappedT &s, jlcxx::ArrayRef<float> m_jl) {
441 0 : std::set<float> m = {m_jl[0], m_jl[1]};
442 0 : s.restrictM(m);
443 : });
444 12 : wrapped.method("diagonalize",
445 : static_cast<void (WrappedT::*)()>(&WrappedT::diagonalize));
446 12 : wrapped.method("getHamiltonian",
447 : static_cast<Eigen::SparseMatrix<Scalar> &(WrappedT::*)()>(
448 : &WrappedT::getHamiltonian));
449 12 : wrapped.method("getOverlap", [](WrappedT &st, StateTwo &s) {
450 0 : Eigen::VectorX<double> overlap = st.getOverlap(s);
451 0 : return jlcxx::get_array_from_evd_t(overlap);
452 : });
453 12 : wrapped.method("getOverlap", [](WrappedT &st, jlcxx::ArrayRef<jl_value_t *> sv) {
454 0 : std::vector<StateTwo> generalizedstates;
455 0 : for (unsigned i = 0; i < sv.size(); i++) {
456 0 : const StateTwo s = *jlcxx::unbox_wrapped_ptr<StateTwo>(sv[i]);
457 0 : generalizedstates.push_back(s);
458 : }
459 0 : const std::vector<StateTwo> &generalizedstates_const = generalizedstates;
460 0 : Eigen::VectorX<double> overlap = st.getOverlap(generalizedstates_const);
461 0 : return jlcxx::get_array_from_evd_t(overlap);
462 : });
463 12 : wrapped.method("getOverlap", [](WrappedT &st, int state_index) {
464 0 : Eigen::VectorX<double> overlap = st.getOverlap(state_index);
465 0 : return jlcxx::get_array_from_evd_t(overlap);
466 : });
467 12 : wrapped.method("getOverlap", [](WrappedT &st, jlcxx::ArrayRef<int> si) {
468 0 : std::vector<size_t> states_indices;
469 0 : for (unsigned i = 0; i < si.size(); i++) {
470 0 : const size_t s = si[i];
471 0 : states_indices.push_back(s);
472 : }
473 0 : const std::vector<size_t> &states_indices_const = states_indices;
474 0 : Eigen::VectorX<double> overlap = st.getOverlap(states_indices_const);
475 0 : return jlcxx::get_array_from_evd_t(overlap);
476 : });
477 12 : wrapped.method("getOverlap",
478 0 : [](WrappedT &st, StateTwo &s, jlcxx::ArrayRef<double> to_z_axis_jl,
479 : jlcxx::ArrayRef<double> to_y_axis_jl) {
480 0 : std::array<double, 3> to_z_axis = {to_z_axis_jl[0], to_z_axis_jl[1],
481 : to_z_axis_jl[2]};
482 0 : std::array<double, 3> to_y_axis = {to_y_axis_jl[0], to_y_axis_jl[1],
483 : to_y_axis_jl[2]};
484 :
485 0 : Eigen::VectorX<double> overlap =
486 : st.getOverlap(s, to_z_axis, to_y_axis);
487 0 : return jlcxx::get_array_from_evd_t(overlap);
488 : });
489 12 : wrapped.method(
490 : "getOverlap",
491 0 : [](WrappedT &st, jlcxx::ArrayRef<jl_value_t *> sv,
492 : jlcxx::ArrayRef<double> to_z_axis_jl, jlcxx::ArrayRef<double> to_y_axis_jl) {
493 0 : std::vector<StateTwo> generalizedstates;
494 0 : for (unsigned i = 0; i < sv.size(); i++) {
495 0 : const StateTwo s = *jlcxx::unbox_wrapped_ptr<StateTwo>(sv[i]);
496 0 : generalizedstates.push_back(s);
497 : }
498 0 : const std::vector<StateTwo> &generalizedstates_const = generalizedstates;
499 :
500 0 : std::array<double, 3> to_z_axis = {to_z_axis_jl[0], to_z_axis_jl[1],
501 : to_z_axis_jl[2]};
502 0 : std::array<double, 3> to_y_axis = {to_y_axis_jl[0], to_y_axis_jl[1],
503 : to_y_axis_jl[2]};
504 :
505 0 : Eigen::VectorX<double> overlap =
506 : st.getOverlap(generalizedstates_const, to_z_axis, to_y_axis);
507 0 : return jlcxx::get_array_from_evd_t(overlap);
508 : });
509 12 : wrapped.method("getOverlap",
510 0 : [](WrappedT &st, int state_index, jlcxx::ArrayRef<double> to_z_axis_jl,
511 : jlcxx::ArrayRef<double> to_y_axis_jl) {
512 0 : std::array<double, 3> to_z_axis = {to_z_axis_jl[0], to_z_axis_jl[1],
513 : to_z_axis_jl[2]};
514 0 : std::array<double, 3> to_y_axis = {to_y_axis_jl[0], to_y_axis_jl[1],
515 : to_y_axis_jl[2]};
516 :
517 0 : Eigen::VectorX<double> overlap =
518 0 : st.getOverlap(state_index, to_z_axis, to_y_axis);
519 0 : return jlcxx::get_array_from_evd_t(overlap);
520 : });
521 12 : wrapped.method("getOverlap",
522 0 : [](WrappedT &st, jlcxx::ArrayRef<int> si,
523 : jlcxx::ArrayRef<double> to_z_axis_jl,
524 : jlcxx::ArrayRef<double> to_y_axis_jl) {
525 0 : std::vector<size_t> states_indices;
526 0 : for (unsigned i = 0; i < si.size(); i++) {
527 0 : const size_t s = si[i];
528 0 : states_indices.push_back(s);
529 : }
530 0 : const std::vector<size_t> &states_indices_const = states_indices;
531 :
532 0 : std::array<double, 3> to_z_axis = {to_z_axis_jl[0], to_z_axis_jl[1],
533 : to_z_axis_jl[2]};
534 0 : std::array<double, 3> to_y_axis = {to_y_axis_jl[0], to_y_axis_jl[1],
535 : to_y_axis_jl[2]};
536 :
537 0 : Eigen::VectorX<double> overlap =
538 : st.getOverlap(states_indices_const, to_z_axis, to_y_axis);
539 0 : return jlcxx::get_array_from_evd_t(overlap);
540 : });
541 12 : wrapped.method("getOverlap",
542 4 : [](WrappedT &st, StateTwo &s, double alpha, double beta, double gamma) {
543 4 : Eigen::VectorX<double> overlap =
544 : st.getOverlap(s, alpha, beta, gamma);
545 8 : return jlcxx::get_array_from_evd_t(overlap);
546 : });
547 12 : wrapped.method(
548 : "getOverlap",
549 0 : [](WrappedT &st, int state_index, double alpha, double beta, double gamma) {
550 0 : Eigen::VectorX<double> overlap = st.getOverlap(state_index, alpha, beta, gamma);
551 0 : return jlcxx::get_array_from_evd_t(overlap);
552 : });
553 12 : wrapped.method("getOverlap",
554 0 : [](WrappedT &st, jlcxx::ArrayRef<jl_value_t *> sv, double alpha,
555 : double beta, double gamma) {
556 0 : std::vector<StateTwo> generalizedstates;
557 0 : for (unsigned i = 0; i < sv.size(); i++) {
558 0 : const StateTwo s = *jlcxx::unbox_wrapped_ptr<StateTwo>(sv[i]);
559 0 : generalizedstates.push_back(s);
560 : }
561 0 : const std::vector<StateTwo> &generalizedstates_const =
562 : generalizedstates;
563 :
564 0 : Eigen::VectorX<double> overlap =
565 : st.getOverlap(generalizedstates_const, alpha, beta, gamma);
566 0 : return jlcxx::get_array_from_evd_t(overlap);
567 : });
568 12 : wrapped.method(
569 : "getOverlap",
570 0 : [](WrappedT &st, jlcxx::ArrayRef<int> si, double alpha, double beta, double gamma) {
571 0 : std::vector<size_t> states_indices;
572 0 : for (unsigned i = 0; i < si.size(); i++) {
573 0 : const size_t s = si[i];
574 0 : states_indices.push_back(s);
575 : }
576 0 : const std::vector<size_t> &states_indices_const = states_indices;
577 :
578 0 : Eigen::VectorX<double> overlap =
579 : st.getOverlap(states_indices_const, alpha, beta, gamma);
580 0 : return jlcxx::get_array_from_evd_t(overlap);
581 : });
582 : /////////////////////////////////////////////////////////////////////////////////////////////
583 12 : wrapped.method("getSpecies", [](WrappedT &s) {
584 0 : std::array<std::string, 2> species_arr = s.getSpecies();
585 0 : return std::make_tuple(species_arr[0], species_arr[1]);
586 : });
587 12 : wrapped.method("getStatesFirst", &WrappedT::getStatesFirst);
588 12 : wrapped.method("getStatesSecond", &WrappedT::getStatesSecond);
589 12 : wrapped.method("enableGreenTensor", &WrappedT::enableGreenTensor);
590 12 : wrapped.method("setSurfaceDistance", &WrappedT::setSurfaceDistance);
591 12 : wrapped.method("setAngle", &WrappedT::setAngle);
592 12 : wrapped.method("setDistance", &WrappedT::setDistance);
593 12 : wrapped.method("setDistanceVector", [](WrappedT &s, jlcxx::ArrayRef<double> dvec) {
594 0 : std::array<double, 3> d = {dvec[0], dvec[1], dvec[2]};
595 0 : s.setDistanceVector(d);
596 : });
597 12 : wrapped.method("setOrder", &WrappedT::setOrder);
598 12 : wrapped.method("setConservedParityUnderPermutation",
599 : &WrappedT::setConservedParityUnderPermutation);
600 12 : wrapped.method("setConservedParityUnderInversion",
601 : &WrappedT::setConservedParityUnderInversion);
602 12 : wrapped.method("setConservedParityUnderReflection",
603 : &WrappedT::setConservedParityUnderReflection);
604 12 : wrapped.method("setConservedMomentaUnderRotation",
605 0 : [](WrappedT &s, jlcxx::ArrayRef<int> momenta_jl) {
606 0 : std::set<int> momenta;
607 0 : for (unsigned i = 0; i < momenta_jl.size(); i++) {
608 0 : momenta.insert(momenta_jl[i]);
609 : }
610 0 : const std::set<int> &momenta_const = momenta;
611 0 : s.setConservedMomentaUnderRotation(momenta_const);
612 : });
613 18 : });
614 :
615 12 : pi.add_type<QuantumDefect>("QuantumDefect")
616 6 : .constructor<std::string const &, int, int, double>()
617 6 : .constructor<std::string const &, int, int, double, std::string const &>()
618 12 : .method("n", [](QuantumDefect &qd) { return qd.n; })
619 12 : .method("l", [](QuantumDefect &qd) { return qd.l; })
620 12 : .method("j", [](QuantumDefect &qd) { return qd.j; })
621 13 : .method("ac", [](QuantumDefect &qd) { return qd.ac; })
622 13 : .method("Z", [](QuantumDefect &qd) { return qd.Z; })
623 13 : .method("a1", [](QuantumDefect &qd) { return qd.a1; })
624 13 : .method("a2", [](QuantumDefect &qd) { return qd.a2; })
625 13 : .method("a3", [](QuantumDefect &qd) { return qd.a3; })
626 13 : .method("a4", [](QuantumDefect &qd) { return qd.a4; })
627 13 : .method("rc", [](QuantumDefect &qd) { return qd.rc; })
628 12 : .method("nstar", [](QuantumDefect &qd) { return qd.nstar; })
629 6 : .method("energy", [](QuantumDefect &qd) { return qd.energy; });
630 6 : }
|