pairinteraction
A Rydberg Interaction Calculator
Basis.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2024 Pairinteraction Developers
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
5
15
16#include <cassert>
17#include <numeric>
18#include <set>
19
20namespace pairinteraction {
21
22template <typename Scalar>
23class BasisAtom;
24
25template <typename Derived>
26void Basis<Derived>::perform_sorter_checks(const std::vector<TransformationType> &labels) const {
27 // Check if the labels are valid sorting labels
28 for (const auto &label : labels) {
29 if (!utils::is_sorting(label)) {
30 throw std::invalid_argument("One of the labels is not a valid sorting label.");
31 }
32 }
33}
34
35template <typename Derived>
37 const std::set<TransformationType> &unique_labels) const {
38 // Check if the states are sorted by the requested labels
39 std::set<TransformationType> unique_labels_present;
40 for (const auto &label : get_transformation().transformation_type) {
41 if (!utils::is_sorting(label) || unique_labels_present.size() >= unique_labels.size()) {
42 break;
43 }
44 unique_labels_present.insert(label);
45 }
46 if (unique_labels != unique_labels_present) {
47 throw std::invalid_argument("The states are not sorted by the requested labels.");
48 }
49
50 // Throw a meaningful error if getting the blocks by energy is requested as this might be a
51 // common mistake
52 if (unique_labels.count(TransformationType::SORT_BY_ENERGY) > 0) {
53 throw std::invalid_argument("States do not store the energy and thus no energy blocks can "
54 "be obtained. Use an energy operator instead.");
55 }
58template <typename Derived>
60 : kets(std::move(kets)), coefficients{{static_cast<Eigen::Index>(this->kets.size()),
61 static_cast<Eigen::Index>(this->kets.size())},
63 if (this->kets.empty()) {
64 throw std::invalid_argument("The basis must contain at least one element.");
65 }
66 state_index_to_quantum_number_f.reserve(this->kets.size());
67 state_index_to_quantum_number_m.reserve(this->kets.size());
68 state_index_to_parity.reserve(this->kets.size());
69 ket_to_ket_index.reserve(this->kets.size());
70 size_t index = 0;
71 for (const auto &ket : this->kets) {
72 state_index_to_quantum_number_f.push_back(ket->get_quantum_number_f());
73 state_index_to_quantum_number_m.push_back(ket->get_quantum_number_m());
74 state_index_to_parity.push_back(ket->get_parity());
75 ket_to_ket_index[ket] = index++;
76 if (ket->get_quantum_number_f() == std::numeric_limits<real_t>::max()) {
77 _has_quantum_number_f = false;
78 }
79 if (ket->get_quantum_number_m() == std::numeric_limits<real_t>::max()) {
80 _has_quantum_number_m = false;
81 }
82 if (ket->get_parity() == Parity::UNKNOWN) {
83 _has_parity = false;
84 }
85 }
86 state_index_to_ket_index.resize(this->kets.size());
87 std::iota(state_index_to_ket_index.begin(), state_index_to_ket_index.end(), 0);
88 ket_index_to_state_index.resize(this->kets.size());
89 std::iota(ket_index_to_state_index.begin(), ket_index_to_state_index.end(), 0);
90 coefficients.matrix.setIdentity();
93template <typename Derived>
95 return _has_quantum_number_f;
96}
97
98template <typename Derived>
100 return _has_quantum_number_m;
103template <typename Derived>
105 return _has_parity;
106}
108template <typename Derived>
109const Derived &Basis<Derived>::derived() const {
110 return static_cast<const Derived &>(*this);
112
113template <typename Derived>
115 return kets;
117
118template <typename Derived>
119const Eigen::SparseMatrix<typename Basis<Derived>::scalar_t, Eigen::RowMajor> &
121 return coefficients.matrix;
122}
123
124template <typename Derived>
125Eigen::SparseMatrix<typename Basis<Derived>::scalar_t, Eigen::RowMajor> &
127 return coefficients.matrix;
128}
129
130template <typename Derived>
131int Basis<Derived>::get_ket_index_from_ket(std::shared_ptr<const ket_t> ket) const {
132 if (ket_to_ket_index.count(ket) == 0) {
133 return -1;
134 }
135 return ket_to_ket_index.at(ket);
136}
137
138template <typename Derived>
140Basis<Derived>::get_amplitudes(std::shared_ptr<const ket_t> ket) const {
141 int ket_index = get_ket_index_from_ket(ket);
142 if (ket_index < 0) {
143 throw std::invalid_argument("The ket does not belong to the basis.");
144 }
145 // The following line is a more efficient alternative to
146 // "get_amplitudes(get_canonical_state_from_ket(ket)).transpose()"
147 return coefficients.matrix.row(ket_index);
148}
149
150template <typename Derived>
151Eigen::SparseMatrix<typename Basis<Derived>::scalar_t, Eigen::RowMajor>
152Basis<Derived>::get_amplitudes(std::shared_ptr<const Derived> other) const {
153 return other->coefficients.matrix.adjoint() * coefficients.matrix;
154}
155
156template <typename Derived>
158Basis<Derived>::get_overlaps(std::shared_ptr<const ket_t> ket) const {
159 return get_amplitudes(ket).cwiseAbs2();
160}
161
162template <typename Derived>
163Eigen::SparseMatrix<typename Basis<Derived>::real_t, Eigen::RowMajor>
164Basis<Derived>::get_overlaps(std::shared_ptr<const Derived> other) const {
165 return get_amplitudes(other).cwiseAbs2();
166}
167
168template <typename Derived>
170 real_t quantum_number_f = state_index_to_quantum_number_f.at(state_index);
171 if (quantum_number_f == std::numeric_limits<real_t>::max()) {
172 throw std::invalid_argument("The state does not have a well-defined quantum number f.");
173 }
174 return quantum_number_f;
175}
176
177template <typename Derived>
179 real_t quantum_number_m = state_index_to_quantum_number_m.at(state_index);
180 if (quantum_number_m == std::numeric_limits<real_t>::max()) {
181 throw std::invalid_argument("The state does not have a well-defined quantum number m.");
182 }
183 return quantum_number_m;
184}
185
186template <typename Derived>
187Parity Basis<Derived>::get_parity(size_t state_index) const {
188 Parity parity = state_index_to_parity.at(state_index);
189 if (parity == Parity::UNKNOWN) {
190 throw std::invalid_argument("The state does not have a well-defined parity.");
191 }
192 return parity;
193}
194
195template <typename Derived>
196std::shared_ptr<const typename Basis<Derived>::ket_t>
197Basis<Derived>::get_corresponding_ket(size_t state_index) const {
198 size_t ket_index = state_index_to_ket_index.at(state_index);
199 if (ket_index == std::numeric_limits<int>::max()) {
200 throw std::invalid_argument("The state does not belong to a ket in a well-defined way.");
201 }
202 return kets[ket_index];
203}
204
205template <typename Derived>
206std::shared_ptr<const typename Basis<Derived>::ket_t>
207Basis<Derived>::get_corresponding_ket(std::shared_ptr<const Derived> /*state*/) const {
208 throw std::runtime_error("Not implemented yet.");
209}
210
211template <typename Derived>
212std::shared_ptr<const Derived> Basis<Derived>::get_state(size_t state_index) const {
213 // Create a copy of the current object
214 auto restricted = std::make_shared<Derived>(derived());
215
216 // Restrict the copy to the state with the largest overlap
217 restricted->coefficients.matrix = restricted->coefficients.matrix.col(state_index);
218
219 std::fill(restricted->ket_index_to_state_index.begin(),
220 restricted->ket_index_to_state_index.end(), std::numeric_limits<int>::max());
221 restricted->ket_index_to_state_index[state_index_to_ket_index[state_index]] = 0;
222
223 restricted->state_index_to_quantum_number_f = {state_index_to_quantum_number_f[state_index]};
224 restricted->state_index_to_quantum_number_m = {state_index_to_quantum_number_m[state_index]};
225 restricted->state_index_to_parity = {state_index_to_parity[state_index]};
226 restricted->state_index_to_ket_index = {state_index_to_ket_index[state_index]};
227
228 restricted->_has_quantum_number_f =
229 restricted->state_index_to_quantum_number_f[0] != std::numeric_limits<real_t>::max();
230 restricted->_has_quantum_number_m =
231 restricted->state_index_to_quantum_number_m[0] != std::numeric_limits<real_t>::max();
232 restricted->_has_parity = restricted->state_index_to_parity[0] != Parity::UNKNOWN;
233
234 return restricted;
235}
236
237template <typename Derived>
238std::shared_ptr<const typename Basis<Derived>::ket_t>
239Basis<Derived>::get_ket(size_t ket_index) const {
240 return kets[ket_index];
241}
242
243template <typename Derived>
244std::shared_ptr<const Derived> Basis<Derived>::get_corresponding_state(size_t ket_index) const {
245 size_t state_index = ket_index_to_state_index.at(ket_index);
246 if (state_index == std::numeric_limits<int>::max()) {
247 throw std::runtime_error("The ket does not belong to a state in a well-defined way.");
248 }
249 return get_state(state_index);
250}
251
252template <typename Derived>
253std::shared_ptr<const Derived>
254Basis<Derived>::get_corresponding_state(std::shared_ptr<const ket_t> ket) const {
255 int ket_index = get_ket_index_from_ket(ket);
256 if (ket_index < 0) {
257 throw std::invalid_argument("The ket does not belong to the basis.");
258 }
259 return get_corresponding_state(ket_index);
260}
261
262template <typename Derived>
263size_t Basis<Derived>::get_corresponding_state_index(size_t ket_index) const {
264 int state_index = ket_index_to_state_index.at(ket_index);
265 if (state_index == std::numeric_limits<int>::max()) {
266 throw std::runtime_error("The ket does not belong to a state in a well-defined way.");
267 }
268 return state_index;
269}
270
271template <typename Derived>
272size_t Basis<Derived>::get_corresponding_state_index(std::shared_ptr<const ket_t> ket) const {
273 int ket_index = get_ket_index_from_ket(ket);
274 if (ket_index < 0) {
275 throw std::invalid_argument("The ket does not belong to the basis.");
276 }
277 return get_corresponding_state_index(ket_index);
278}
279
280template <typename Derived>
281size_t Basis<Derived>::get_corresponding_ket_index(size_t state_index) const {
282 int ket_index = state_index_to_ket_index.at(state_index);
283 if (ket_index == std::numeric_limits<int>::max()) {
284 throw std::runtime_error("The state does not belong to a ket in a well-defined way.");
285 }
286 return ket_index;
287}
288
289template <typename Derived>
290size_t Basis<Derived>::get_corresponding_ket_index(std::shared_ptr<const Derived> /*state*/) const {
291 throw std::runtime_error("Not implemented yet.");
292}
293
294template <typename Derived>
295std::shared_ptr<const Derived>
297 // Create a copy of the current object
298 auto created = std::make_shared<Derived>(derived());
299
300 // Fill the copy with the state corresponding to the ket index
301 created->coefficients.matrix =
302 Eigen::SparseMatrix<scalar_t, Eigen::RowMajor>(coefficients.matrix.rows(), 1);
303 created->coefficients.matrix.coeffRef(ket_index, 0) = 1;
304 created->coefficients.matrix.makeCompressed();
305
306 std::fill(created->ket_index_to_state_index.begin(), created->ket_index_to_state_index.end(),
307 std::numeric_limits<int>::max());
308 created->ket_index_to_state_index[ket_index] = 0;
309
310 created->state_index_to_quantum_number_f = {kets[ket_index]->get_quantum_number_f()};
311 created->state_index_to_quantum_number_m = {kets[ket_index]->get_quantum_number_m()};
312 created->state_index_to_parity = {kets[ket_index]->get_parity()};
313 created->state_index_to_ket_index = {ket_index};
314
315 created->_has_quantum_number_f =
316 created->state_index_to_quantum_number_f[0] != std::numeric_limits<real_t>::max();
317 created->_has_quantum_number_m =
318 created->state_index_to_quantum_number_m[0] != std::numeric_limits<real_t>::max();
319 created->_has_parity = created->state_index_to_parity[0] != Parity::UNKNOWN;
320
321 return created;
322}
323
324template <typename Derived>
325std::shared_ptr<const Derived>
326Basis<Derived>::get_canonical_state_from_ket(std::shared_ptr<const ket_t> ket) const {
327 int ket_index = get_ket_index_from_ket(ket);
328 if (ket_index < 0) {
329 throw std::invalid_argument("The ket does not belong to the basis.");
330 }
331 return get_canonical_state_from_ket(ket_index);
332}
333
334template <typename Derived>
336 return kets.begin();
337}
338
339template <typename Derived>
341 return kets.end();
342}
343
344template <typename Derived>
345Basis<Derived>::Iterator::Iterator(typename ketvec_t::const_iterator it) : it{std::move(it)} {}
346
347template <typename Derived>
349 return other.it != it;
350}
351
352template <typename Derived>
353std::shared_ptr<const typename Basis<Derived>::ket_t> Basis<Derived>::Iterator::operator*() const {
354 return *it;
355}
356
357template <typename Derived>
359 ++it;
360 return *this;
361}
362
363template <typename Derived>
365 return coefficients.matrix.cols();
366}
367
368template <typename Derived>
370 return coefficients.matrix.rows();
371}
372
373template <typename Derived>
376 return coefficients;
377}
378
379template <typename Derived>
382 Transformation<scalar_t> transformation{{static_cast<Eigen::Index>(coefficients.matrix.rows()),
383 static_cast<Eigen::Index>(coefficients.matrix.rows())},
385
386 std::vector<Eigen::Triplet<scalar_t>> entries;
387
388 for (size_t idx_initial = 0; idx_initial < kets.size(); ++idx_initial) {
389 real_t f = kets[idx_initial]->get_quantum_number_f();
390 real_t m_initial = kets[idx_initial]->get_quantum_number_m();
391
392 assert(2 * f == std::floor(2 * f) && f >= 0);
393 assert(2 * m_initial == std::floor(2 * m_initial) && m_initial >= -f && m_initial <= f);
394
395 for (real_t m_final = -f; m_final <= f; // NOSONAR m_final is precisely representable
396 ++m_final) {
397 auto val = wigner::wigner_uppercase_d_matrix<scalar_t>(f, m_initial, m_final, alpha,
398 beta, gamma);
399 size_t idx_final = get_ket_index_from_ket(
400 kets[idx_initial]->get_ket_for_different_quantum_number_m(m_final));
401 entries.emplace_back(idx_final, idx_initial, val);
402 }
403 }
404
405 transformation.matrix.setFromTriplets(entries.begin(), entries.end());
406 transformation.matrix.makeCompressed();
407
408 return transformation;
409}
410
411template <typename Derived>
412Sorting Basis<Derived>::get_sorter(const std::vector<TransformationType> &labels) const {
413 perform_sorter_checks(labels);
414
415 // Throw a meaningful error if sorting by energy is requested as this might be a common mistake
416 if (std::find(labels.begin(), labels.end(), TransformationType::SORT_BY_ENERGY) !=
417 labels.end()) {
418 throw std::invalid_argument("States do not store the energy and thus can not be sorted by "
419 "the energy. Use an energy operator instead.");
420 }
421
422 // Initialize transformation
423 Sorting transformation;
424 transformation.matrix.resize(coefficients.matrix.cols());
425 transformation.matrix.setIdentity();
426
427 // Get the sorter
428 get_sorter_without_checks(labels, transformation);
429
430 // Check if all labels have been used for sorting
431 if (labels != transformation.transformation_type) {
432 throw std::invalid_argument("The states could not be sorted by all the requested labels.");
433 }
434
435 return transformation;
436}
437
438template <typename Derived>
439std::vector<IndicesOfBlock>
440Basis<Derived>::get_indices_of_blocks(const std::vector<TransformationType> &labels) const {
441 perform_sorter_checks(labels);
442
443 std::set<TransformationType> unique_labels(labels.begin(), labels.end());
444 perform_blocks_checks(unique_labels);
445
446 // Get the blocks
447 IndicesOfBlocksCreator blocks_creator({0, static_cast<size_t>(coefficients.matrix.cols())});
448 get_indices_of_blocks_without_checks(unique_labels, blocks_creator);
449
450 return blocks_creator.create();
451}
452
453template <typename Derived>
454void Basis<Derived>::get_sorter_without_checks(const std::vector<TransformationType> &labels,
455 Sorting &transformation) const {
456 constexpr real_t numerical_precision = 100 * std::numeric_limits<real_t>::epsilon();
457
458 int *perm_begin = transformation.matrix.indices().data();
459 int *perm_end = perm_begin + coefficients.matrix.cols();
460 const int *perm_back = perm_end - 1;
461
462 // Sort the vector based on the requested labels
463 std::stable_sort(perm_begin, perm_end, [&](int a, int b) {
464 for (const auto &label : labels) {
465 switch (label) {
467 if (state_index_to_parity[a] != state_index_to_parity[b]) {
468 return state_index_to_parity[a] < state_index_to_parity[b];
469 }
470 break;
472 if (std::abs(state_index_to_quantum_number_m[a] -
473 state_index_to_quantum_number_m[b]) > numerical_precision) {
474 return state_index_to_quantum_number_m[a] < state_index_to_quantum_number_m[b];
475 }
476 break;
478 if (std::abs(state_index_to_quantum_number_f[a] -
479 state_index_to_quantum_number_f[b]) > numerical_precision) {
480 return state_index_to_quantum_number_f[a] < state_index_to_quantum_number_f[b];
481 }
482 break;
484 if (state_index_to_ket_index[a] != state_index_to_ket_index[b]) {
485 return state_index_to_ket_index[a] < state_index_to_ket_index[b];
486 }
487 break;
488 default:
489 std::abort(); // Can't happen because of previous checks
490 }
491 }
492 return false; // Elements are equal
493 });
494
495 // Check for invalid values and add transformation types
496 for (const auto &label : labels) {
497 switch (label) {
499 if (state_index_to_parity[*perm_back] == Parity::UNKNOWN) {
500 throw std::invalid_argument(
501 "States cannot be labeled and thus not sorted by the parity.");
502 }
504 break;
506 if (state_index_to_quantum_number_m[*perm_back] == std::numeric_limits<real_t>::max()) {
507 throw std::invalid_argument(
508 "States cannot be labeled and thus not sorted by the quantum number m.");
509 }
510 transformation.transformation_type.push_back(
512 break;
514 if (state_index_to_quantum_number_f[*perm_back] == std::numeric_limits<real_t>::max()) {
515 throw std::invalid_argument(
516 "States cannot be labeled and thus not sorted by the quantum number f.");
517 }
518 transformation.transformation_type.push_back(
520 break;
522 if (state_index_to_ket_index[*perm_back] == std::numeric_limits<int>::max()) {
523 throw std::invalid_argument(
524 "States cannot be labeled and thus not sorted by kets.");
525 }
527 break;
528 default:
529 std::abort(); // Can't happen because of previous checks
530 }
531 }
532}
533
534template <typename Derived>
536 const std::set<TransformationType> &unique_labels,
537 IndicesOfBlocksCreator &blocks_creator) const {
538 constexpr real_t numerical_precision = 100 * std::numeric_limits<real_t>::epsilon();
539
540 auto last_quantum_number_f = state_index_to_quantum_number_f[0];
541 auto last_quantum_number_m = state_index_to_quantum_number_m[0];
542 auto last_parity = state_index_to_parity[0];
543 auto last_ket = state_index_to_ket_index[0];
544
545 for (int i = 0; i < coefficients.matrix.cols(); ++i) {
546 for (auto label : unique_labels) {
548 if (std::abs(state_index_to_quantum_number_f[i] - last_quantum_number_f) >
549 numerical_precision) {
550 blocks_creator.add(i);
551 break;
552 }
554 if (std::abs(state_index_to_quantum_number_m[i] - last_quantum_number_m) >
555 numerical_precision) {
556 blocks_creator.add(i);
557 break;
558 }
559 } else if (label == TransformationType::SORT_BY_PARITY) {
560 if (state_index_to_parity[i] != last_parity) {
561 blocks_creator.add(i);
562 break;
563 }
564 } else if (label == TransformationType::SORT_BY_KET) {
565 if (state_index_to_ket_index[i] != last_ket) {
566 blocks_creator.add(i);
567 break;
568 }
569 }
570 }
571 last_quantum_number_f = state_index_to_quantum_number_f[i];
572 last_quantum_number_m = state_index_to_quantum_number_m[i];
573 last_parity = state_index_to_parity[i];
574 last_ket = state_index_to_ket_index[i];
575 }
576}
577
578template <typename Derived>
579std::shared_ptr<const Derived> Basis<Derived>::transformed(const Sorting &transformation) const {
580 // Create a copy of the current object
581 auto transformed = std::make_shared<Derived>(derived());
582
583 if (coefficients.matrix.cols() == 0) {
584 return transformed;
585 }
586
587 // Apply the transformation
588 transformed->coefficients.matrix = coefficients.matrix * transformation.matrix;
589 transformed->coefficients.transformation_type = transformation.transformation_type;
590
591 transformed->state_index_to_quantum_number_f.resize(transformation.matrix.size());
592 transformed->state_index_to_quantum_number_m.resize(transformation.matrix.size());
593 transformed->state_index_to_parity.resize(transformation.matrix.size());
594 transformed->state_index_to_ket_index.resize(transformation.matrix.size());
595
596 for (int i = 0; i < transformation.matrix.size(); ++i) {
597 transformed->state_index_to_quantum_number_f[i] =
598 state_index_to_quantum_number_f[transformation.matrix.indices()[i]];
599 transformed->state_index_to_quantum_number_m[i] =
600 state_index_to_quantum_number_m[transformation.matrix.indices()[i]];
601 transformed->state_index_to_parity[i] =
602 state_index_to_parity[transformation.matrix.indices()[i]];
603 transformed->state_index_to_ket_index[i] =
604 state_index_to_ket_index[transformation.matrix.indices()[i]];
605 transformed->ket_index_to_state_index
606 [state_index_to_ket_index[transformation.matrix.indices()[i]]] = i;
607 }
608
609 return transformed;
610}
611
612template <typename Derived>
613std::shared_ptr<const Derived>
615 // TODO why is "numerical_precision = 100 * std::sqrt(coefficients.matrix.rows()) *
616 // std::numeric_limits<real_t>::epsilon()" too small for figuring out whether m is conserved?
617 real_t numerical_precision = 0.001;
618
619 // If the transformation is a rotation, it should be a rotation and nothing else
620 bool is_rotation = false;
621 for (auto t : transformation.transformation_type) {
623 is_rotation = true;
624 break;
625 }
626 }
627 if (is_rotation && transformation.transformation_type.size() != 1) {
628 throw std::invalid_argument("A rotation can not be combined with other transformations.");
629 }
630
631 // To apply a rotation, the object must only be sorted but other transformations are not allowed
632 if (is_rotation) {
633 for (auto t : coefficients.transformation_type) {
634 if (!utils::is_sorting(t)) {
635 throw std::runtime_error(
636 "If the object was transformed by a different transformation "
637 "than sorting, it can not be rotated.");
638 }
639 }
640 }
641
642 // Create a copy of the current object
643 auto transformed = std::make_shared<Derived>(derived());
644
645 if (coefficients.matrix.cols() == 0) {
646 return transformed;
647 }
648
649 // Apply the transformation
650 // If a quantum number turns out to be conserved by the transformation, it will be
651 // rounded to the nearest half integer to avoid loss of numerical_precision.
652 transformed->coefficients.matrix = coefficients.matrix * transformation.matrix;
653 transformed->coefficients.transformation_type = transformation.transformation_type;
654
655 Eigen::SparseMatrix<real_t> probs = transformation.matrix.cwiseAbs2().transpose();
656
657 {
658 auto map = Eigen::Map<const Eigen::VectorX<real_t>>(state_index_to_quantum_number_f.data(),
659 state_index_to_quantum_number_f.size());
660 Eigen::VectorX<real_t> val = probs * map;
661 Eigen::VectorX<real_t> sq = probs * map.cwiseAbs2();
662 Eigen::VectorX<real_t> diff = (val.cwiseAbs2() - sq).cwiseAbs();
663 transformed->state_index_to_quantum_number_f.resize(probs.rows());
664
665 for (size_t i = 0; i < transformed->state_index_to_quantum_number_f.size(); ++i) {
666 if (diff[i] < numerical_precision) {
667 transformed->state_index_to_quantum_number_f[i] = std::round(val[i] * 2) / 2;
668 } else {
669 transformed->state_index_to_quantum_number_f[i] =
670 std::numeric_limits<real_t>::max();
671 transformed->_has_quantum_number_f = false;
672 }
673 }
674 }
675
676 {
677 auto map = Eigen::Map<const Eigen::VectorX<real_t>>(state_index_to_quantum_number_m.data(),
678 state_index_to_quantum_number_m.size());
679 Eigen::VectorX<real_t> val = probs * map;
680 Eigen::VectorX<real_t> sq = probs * map.cwiseAbs2();
681 Eigen::VectorX<real_t> diff = (val.cwiseAbs2() - sq).cwiseAbs();
682 transformed->state_index_to_quantum_number_m.resize(probs.rows());
683
684 for (size_t i = 0; i < transformed->state_index_to_quantum_number_m.size(); ++i) {
685 if (diff[i] < numerical_precision) {
686 transformed->state_index_to_quantum_number_m[i] = std::round(val[i] * 2) / 2;
687 } else {
688 transformed->state_index_to_quantum_number_m[i] =
689 std::numeric_limits<real_t>::max();
690 transformed->_has_quantum_number_m = false;
691 }
692 }
693 }
694
695 {
696 using utype = std::underlying_type<Parity>::type;
697 Eigen::VectorX<real_t> map(state_index_to_parity.size());
698 for (size_t i = 0; i < state_index_to_parity.size(); ++i) {
699 map[i] = static_cast<utype>(state_index_to_parity[i]);
700 }
701 Eigen::VectorX<real_t> val = probs * map;
702 Eigen::VectorX<real_t> sq = probs * map.cwiseAbs2();
703 Eigen::VectorX<real_t> diff = (val.cwiseAbs2() - sq).cwiseAbs();
704 transformed->state_index_to_parity.resize(probs.rows());
705
706 for (size_t i = 0; i < transformed->state_index_to_parity.size(); ++i) {
707 if (diff[i] < numerical_precision) {
708 transformed->state_index_to_parity[i] = static_cast<Parity>(std::lround(val[i]));
709 } else {
710 transformed->state_index_to_parity[i] = Parity::UNKNOWN;
711 transformed->_has_parity = false;
712 }
713 }
714 }
715
716 {
717 // In the following, we obtain a bijective map between state index and ket index.
718
719 // Find the maximum value in each row and column
720 std::vector<real_t> max_in_row(transformed->coefficients.matrix.rows(), 0);
721 std::vector<real_t> max_in_col(transformed->coefficients.matrix.cols(), 0);
722 for (int row = 0; row < transformed->coefficients.matrix.outerSize(); ++row) {
723 for (typename Eigen::SparseMatrix<scalar_t, Eigen::RowMajor>::InnerIterator it(
724 transformed->coefficients.matrix, row);
725 it; ++it) {
726 real_t val = std::pow(std::abs(it.value()), 2);
727 max_in_row[row] = std::max(max_in_row[row], val);
728 max_in_col[it.col()] = std::max(max_in_col[it.col()], val);
729 }
730 }
731
732 // Use the maximum values to define a cost for a sub-optimal mapping
733 std::vector<real_t> costs;
734 std::vector<std::pair<int, int>> mappings;
735 costs.reserve(transformed->coefficients.matrix.nonZeros());
736 mappings.reserve(transformed->coefficients.matrix.nonZeros());
737 for (int row = 0; row < transformed->coefficients.matrix.outerSize(); ++row) {
738 for (typename Eigen::SparseMatrix<scalar_t, Eigen::RowMajor>::InnerIterator it(
739 transformed->coefficients.matrix, row);
740 it; ++it) {
741 real_t val = std::pow(std::abs(it.value()), 2);
742 real_t cost = max_in_row[row] + max_in_col[it.col()] - 2 * val;
743 costs.push_back(cost);
744 mappings.push_back({row, it.col()});
745 }
746 }
747
748 // Obtain from the costs in which order the mappings should be considered
749 std::vector<size_t> order(costs.size());
750 std::iota(order.begin(), order.end(), 0);
751 std::sort(order.begin(), order.end(),
752 [&](size_t a, size_t b) { return costs[a] < costs[b]; });
753
754 // Fill ket_index_to_state_index with invalid values as there can be more kets than states
755 std::fill(transformed->ket_index_to_state_index.begin(),
756 transformed->ket_index_to_state_index.end(), std::numeric_limits<int>::max());
757
758 // Generate the bijective map
759 std::vector<bool> row_used(transformed->coefficients.matrix.rows(), false);
760 std::vector<bool> col_used(transformed->coefficients.matrix.cols(), false);
761 int num_used = 0;
762 for (size_t idx : order) {
763 int row = mappings[idx].first; // corresponds to the ket index
764 int col = mappings[idx].second; // corresponds to the state index
765 if (!row_used[row] && !col_used[col]) {
766 row_used[row] = true;
767 col_used[col] = true;
768 num_used++;
769 transformed->state_index_to_ket_index[col] = row;
770 transformed->ket_index_to_state_index[row] = col;
771 }
772 if (num_used == transformed->coefficients.matrix.cols()) {
773 break;
774 }
775 }
776 assert(num_used == transformed->coefficients.matrix.cols());
777 }
778
779 return transformed;
780}
781
782template <typename Derived>
783size_t Basis<Derived>::hash::operator()(const std::shared_ptr<const ket_t> &k) const {
784 return typename ket_t::hash()(*k);
785}
786
787template <typename Derived>
788bool Basis<Derived>::equal_to::operator()(const std::shared_ptr<const ket_t> &lhs,
789 const std::shared_ptr<const ket_t> &rhs) const {
790 return *lhs == *rhs;
791}
792
793// Explicit instantiations
794template class Basis<BasisAtom<double>>;
795template class Basis<BasisAtom<std::complex<double>>>;
796template class Basis<BasisPair<double>>;
797template class Basis<BasisPair<std::complex<double>>>;
798} // namespace pairinteraction
std::shared_ptr< const ket_t > operator*() const
Definition: Basis.cpp:353
bool operator!=(const Iterator &other) const
Definition: Basis.cpp:348
Base class for a basis.
Definition: Basis.hpp:41
const Transformation< scalar_t > & get_transformation() const override
Definition: Basis.cpp:375
Iterator begin() const
Definition: Basis.cpp:335
void get_indices_of_blocks_without_checks(const std::set< TransformationType > &unique_labels, IndicesOfBlocksCreator &blocks) const
Definition: Basis.cpp:535
void get_sorter_without_checks(const std::vector< TransformationType > &labels, Sorting &transformation) const
Definition: Basis.cpp:454
int get_ket_index_from_ket(std::shared_ptr< const ket_t > ket) const
Definition: Basis.cpp:131
typename traits::CrtpTraits< Derived >::ketvec_t ketvec_t
Definition: Basis.hpp:46
Iterator end() const
Definition: Basis.cpp:340
size_t get_number_of_states() const
Definition: Basis.cpp:364
Transformation< scalar_t > get_rotator(real_t alpha, real_t beta, real_t gamma) const override
Definition: Basis.cpp:381
void perform_sorter_checks(const std::vector< TransformationType > &labels) const
Definition: Basis.cpp:26
void perform_blocks_checks(const std::set< TransformationType > &unique_labels) const
Definition: Basis.cpp:36
std::vector< IndicesOfBlock > get_indices_of_blocks(const std::vector< TransformationType > &labels) const override
Definition: Basis.cpp:440
typename traits::CrtpTraits< Derived >::real_t real_t
Definition: Basis.hpp:44
std::shared_ptr< const Derived > transformed(const Transformation< scalar_t > &transformation) const
Definition: Basis.cpp:614
Sorting get_sorter(const std::vector< TransformationType > &labels) const override
Definition: Basis.cpp:412
size_t get_number_of_kets() const
Definition: Basis.cpp:369
Matrix< Type, Dynamic, 1 > VectorX
bool is_sorting(TransformationType label)
Eigen::PermutationMatrix< Eigen::Dynamic, Eigen::Dynamic > matrix
std::vector< TransformationType > transformation_type
Eigen::SparseMatrix< Scalar, Eigen::RowMajor > matrix
std::vector< TransformationType > transformation_type