Line data Source code
1 : // SPDX-FileCopyrightText: 2024 PairInteraction Developers
2 : // SPDX-License-Identifier: LGPL-3.0-or-later
3 :
4 : #include "pairinteraction/system/System.hpp"
5 :
6 : #include "pairinteraction/basis/BasisAtom.hpp"
7 : #include "pairinteraction/basis/BasisPair.hpp"
8 : #include "pairinteraction/enums/TransformationType.hpp"
9 : #include "pairinteraction/interfaces/DiagonalizerInterface.hpp"
10 : #include "pairinteraction/system/SystemAtom.hpp"
11 : #include "pairinteraction/system/SystemPair.hpp"
12 : #include "pairinteraction/utils/TaskControl.hpp"
13 : #include "pairinteraction/utils/eigen_assertion.hpp"
14 : #include "pairinteraction/utils/eigen_compat.hpp"
15 :
16 : #include <Eigen/SparseCore>
17 : #include <algorithm>
18 : #include <complex>
19 : #include <limits>
20 : #include <memory>
21 : #include <numeric>
22 : #include <oneapi/tbb.h>
23 : #include <optional>
24 : #include <spdlog/spdlog.h>
25 :
26 : namespace pairinteraction {
27 : template <typename Derived>
28 845 : System<Derived>::System(std::shared_ptr<const basis_t> basis)
29 845 : : basis(std::move(basis)),
30 845 : matrix(static_cast<Eigen::Index>(this->basis->get_number_of_states()),
31 2535 : static_cast<Eigen::Index>(this->basis->get_number_of_states())) {}
32 :
33 : template <typename Derived>
34 581 : std::shared_ptr<const typename System<Derived>::basis_t> System<Derived>::get_basis() const {
35 581 : if (hamiltonian_requires_construction) {
36 17 : construct_hamiltonian();
37 17 : hamiltonian_requires_construction = false;
38 : }
39 581 : return basis;
40 : }
41 :
42 : template <typename Derived>
43 261 : std::shared_ptr<const typename System<Derived>::basis_t> System<Derived>::get_eigenbasis() const {
44 261 : if (hamiltonian_requires_construction) {
45 0 : construct_hamiltonian();
46 0 : hamiltonian_requires_construction = false;
47 : }
48 261 : if (!this->is_diagonal()) {
49 0 : throw std::runtime_error("The Hamiltonian has not been diagonalized yet.");
50 : }
51 261 : return basis;
52 : }
53 :
54 : template <typename Derived>
55 1190 : Eigen::VectorX<typename System<Derived>::real_t> System<Derived>::get_eigenenergies() const {
56 1190 : if (hamiltonian_requires_construction) {
57 36 : construct_hamiltonian();
58 36 : hamiltonian_requires_construction = false;
59 : }
60 1190 : if (!this->is_diagonal()) {
61 0 : throw std::runtime_error("The Hamiltonian has not been diagonalized yet.");
62 : }
63 1864 : return matrix.diagonal().real();
64 : }
65 :
66 : template <typename Derived>
67 : const Eigen::SparseMatrix<typename System<Derived>::scalar_t, Eigen::RowMajor> &
68 154 : System<Derived>::get_matrix() const {
69 154 : if (hamiltonian_requires_construction) {
70 113 : construct_hamiltonian();
71 113 : hamiltonian_requires_construction = false;
72 : }
73 154 : return matrix;
74 : }
75 :
76 : template <typename Derived>
77 : const Transformation<typename System<Derived>::scalar_t> &
78 0 : System<Derived>::get_transformation() const {
79 0 : if (hamiltonian_requires_construction) {
80 0 : construct_hamiltonian();
81 0 : hamiltonian_requires_construction = false;
82 : }
83 0 : return basis->get_transformation();
84 : }
85 :
86 : template <typename Derived>
87 : Transformation<typename System<Derived>::scalar_t>
88 0 : System<Derived>::get_rotator(real_t alpha, real_t beta, real_t gamma) const {
89 0 : if (hamiltonian_requires_construction) {
90 0 : construct_hamiltonian();
91 0 : hamiltonian_requires_construction = false;
92 : }
93 0 : return basis->get_rotator(alpha, beta, gamma);
94 : }
95 :
96 : template <typename Derived>
97 1615 : Sorting System<Derived>::get_sorter(const std::vector<TransformationType> &labels) const {
98 1615 : if (hamiltonian_requires_construction) {
99 0 : construct_hamiltonian();
100 0 : hamiltonian_requires_construction = false;
101 : }
102 :
103 1615 : basis->perform_sorter_checks(labels);
104 :
105 1615 : auto it = std::find(labels.begin(), labels.end(), TransformationType::SORT_BY_ENERGY);
106 1615 : std::vector<TransformationType> before_energy(labels.begin(), it);
107 1615 : bool contains_energy = (it != labels.end());
108 1615 : std::vector<TransformationType> after_energy(contains_energy ? it + 1 : labels.end(),
109 : labels.end());
110 :
111 1615 : Sorting transformation;
112 1615 : transformation.matrix.resize(matrix.rows());
113 1615 : transformation.matrix.setIdentity();
114 :
115 1615 : if (!before_energy.empty()) {
116 568 : basis->get_sorter_without_checks(before_energy, transformation);
117 : }
118 :
119 1615 : if (contains_energy) {
120 1047 : std::vector<real_t> energies_of_states;
121 1047 : energies_of_states.reserve(matrix.rows());
122 120813 : for (int i = 0; i < matrix.rows(); ++i) {
123 119766 : energies_of_states.push_back(std::real(matrix.coeff(i, i)));
124 : }
125 :
126 1047 : std::stable_sort(
127 1047 : transformation.matrix.indices().data(),
128 1047 : transformation.matrix.indices().data() + transformation.matrix.indices().size(),
129 555822 : [&](int i, int j) { return energies_of_states[i] < energies_of_states[j]; });
130 :
131 1047 : transformation.transformation_type.push_back(TransformationType::SORT_BY_ENERGY);
132 1047 : }
133 :
134 1615 : if (!after_energy.empty()) {
135 0 : basis->get_sorter_without_checks(after_energy, transformation);
136 : }
137 :
138 1615 : if (labels != transformation.transformation_type) {
139 0 : throw std::invalid_argument("The states could not be sorted by all the requested labels.");
140 : }
141 :
142 3230 : return transformation;
143 1615 : }
144 :
145 : template <typename Derived>
146 : std::vector<IndicesOfBlock>
147 621 : System<Derived>::get_indices_of_blocks(const std::vector<TransformationType> &labels) const {
148 621 : if (hamiltonian_requires_construction) {
149 0 : construct_hamiltonian();
150 0 : hamiltonian_requires_construction = false;
151 : }
152 :
153 621 : basis->perform_sorter_checks(labels);
154 :
155 621 : std::set<TransformationType> unique_labels(labels.begin(), labels.end());
156 621 : basis->perform_blocks_checks(unique_labels);
157 :
158 621 : auto it = unique_labels.find(TransformationType::SORT_BY_ENERGY);
159 621 : bool contains_energy = (it != unique_labels.end());
160 621 : if (contains_energy) {
161 0 : unique_labels.erase(it);
162 : }
163 :
164 621 : IndicesOfBlocksCreator blocks_creator({0, static_cast<size_t>(matrix.rows())});
165 :
166 621 : if (!unique_labels.empty()) {
167 568 : basis->get_indices_of_blocks_without_checks(unique_labels, blocks_creator);
168 : }
169 :
170 621 : if (contains_energy && matrix.rows() > 0) {
171 0 : scalar_t last_energy = std::real(matrix.coeff(0, 0));
172 0 : for (int i = 0; i < matrix.rows(); ++i) {
173 0 : if (std::real(matrix.coeff(i, i)) != last_energy) {
174 0 : blocks_creator.add(i);
175 0 : last_energy = std::real(matrix.coeff(i, i));
176 : }
177 : }
178 : }
179 :
180 1242 : return blocks_creator.create();
181 621 : }
182 :
183 : template <typename Derived>
184 0 : System<Derived> &System<Derived>::transform(const Transformation<scalar_t> &transformation) {
185 0 : if (hamiltonian_requires_construction) {
186 0 : construct_hamiltonian();
187 0 : hamiltonian_requires_construction = false;
188 : }
189 :
190 0 : set_task_status("Applying transformation...");
191 0 : if (matrix.cols() != 0) {
192 0 : matrix = transformation.matrix.adjoint() * matrix * transformation.matrix;
193 0 : basis = basis->transformed(transformation);
194 : }
195 :
196 0 : hamiltonian_is_diagonal = false;
197 :
198 : // A transformed system might have lost its block-diagonalizability if the
199 : // transformation was not a sorting
200 0 : blockdiagonalizing_labels.clear();
201 :
202 0 : return *this;
203 : }
204 :
205 : template <typename Derived>
206 1047 : System<Derived> &System<Derived>::transform(const Sorting &transformation) {
207 1047 : if (hamiltonian_requires_construction) {
208 0 : construct_hamiltonian();
209 0 : hamiltonian_requires_construction = false;
210 : }
211 :
212 1047 : set_task_status("Applying transformation...");
213 1047 : if (matrix.cols() != 0) {
214 1047 : matrix = matrix.twistedBy(transformation.matrix.inverse());
215 1047 : basis = basis->transformed(transformation);
216 : }
217 :
218 1047 : return *this;
219 : }
220 :
221 : template <typename Derived>
222 652 : System<Derived> &System<Derived>::diagonalize(const DiagonalizerInterface<scalar_t> &diagonalizer,
223 : std::optional<real_t> min_eigenenergy,
224 : std::optional<real_t> max_eigenenergy, double rtol) {
225 652 : set_task_status("Preparing Hamiltonian...");
226 :
227 652 : if (hamiltonian_requires_construction) {
228 638 : construct_hamiltonian();
229 638 : hamiltonian_requires_construction = false;
230 : }
231 :
232 652 : if (this->is_diagonal()) {
233 31 : return *this;
234 : }
235 :
236 621 : Eigen::SparseMatrix<scalar_t, Eigen::RowMajor> eigenvectors;
237 621 : Eigen::SparseMatrix<scalar_t, Eigen::RowMajor> eigenenergies;
238 :
239 : // Sort the Hamiltonian according to the block structure
240 621 : if (!blockdiagonalizing_labels.empty()) {
241 568 : auto sorter = get_sorter(blockdiagonalizing_labels);
242 568 : matrix = matrix.twistedBy(sorter.matrix.inverse());
243 568 : basis = basis->transformed(sorter);
244 568 : }
245 :
246 : // Get the indices of the blocks
247 621 : auto blocks = get_indices_of_blocks(blockdiagonalizing_labels);
248 :
249 621 : assert((blockdiagonalizing_labels.empty() && blocks.size() == 1) ||
250 : !blockdiagonalizing_labels.empty());
251 :
252 621 : SPDLOG_DEBUG("Diagonalizing the Hamiltonian with {} blocks.", blocks.size());
253 :
254 : // Diagonalize the blocks in parallel
255 621 : std::vector<Eigen::VectorX<real_t>> eigenenergies_blocks(blocks.size());
256 621 : std::vector<Eigen::SparseMatrix<scalar_t, Eigen::RowMajor>> eigenvectors_blocks(blocks.size());
257 621 : oneapi::tbb::parallel_for(
258 1938 : oneapi::tbb::blocked_range<size_t>(0, blocks.size()), [&](const auto &range) {
259 2634 : for (size_t idx = range.begin(); idx != range.end(); ++idx) {
260 1317 : set_task_status("Diagonalizing Hamiltonian blocks...");
261 4316 : auto eigensys = min_eigenenergy.has_value() || max_eigenenergy.has_value()
262 2047 : ? diagonalizer.eigh(matrix.block(blocks[idx].start, blocks[idx].start,
263 730 : blocks[idx].size(), blocks[idx].size()),
264 : min_eigenenergy, max_eigenenergy, rtol)
265 3221 : : diagonalizer.eigh(matrix.block(blocks[idx].start, blocks[idx].start,
266 1904 : blocks[idx].size(), blocks[idx].size()),
267 : rtol);
268 1317 : eigenvectors_blocks[idx] = eigensys.eigenvectors;
269 1317 : eigenenergies_blocks[idx] = eigensys.eigenvalues;
270 : }
271 : });
272 :
273 : // Get the number of non-zeros per row of the combined eigenvector matrix
274 621 : std::vector<Eigen::Index> non_zeros_per_inner_index;
275 621 : non_zeros_per_inner_index.reserve(matrix.rows());
276 621 : Eigen::Index num_rows = 0;
277 621 : Eigen::Index num_cols = 0;
278 1938 : for (const auto &matrix : eigenvectors_blocks) {
279 1317 : set_task_status("Collecting eigenvectors...");
280 36629 : for (int i = 0; i < matrix.outerSize(); ++i) {
281 70624 : non_zeros_per_inner_index.push_back(matrix.outerIndexPtr()[i + 1] -
282 35312 : matrix.outerIndexPtr()[i]);
283 : }
284 1317 : num_rows += matrix.rows();
285 1317 : num_cols += matrix.cols();
286 : }
287 :
288 621 : assert(static_cast<size_t>(num_rows) == basis->get_number_of_kets());
289 621 : assert(static_cast<size_t>(num_cols) <= basis->get_number_of_states());
290 :
291 621 : eigenvectors.resize(num_rows, num_cols);
292 621 : eigenenergies.resize(num_cols, num_cols);
293 :
294 621 : if (num_cols > 0) {
295 : // Get the combined eigenvector matrix (in case of an restricted energy range, it is not
296 : // square)
297 619 : eigenvectors.reserve(non_zeros_per_inner_index);
298 619 : Eigen::Index offset_rows = 0;
299 619 : Eigen::Index offset_cols = 0;
300 1934 : for (const auto &matrix : eigenvectors_blocks) {
301 1315 : set_task_status("Combining eigenvectors...");
302 36547 : for (Eigen::Index i = 0; i < matrix.outerSize(); ++i) {
303 35232 : for (typename Eigen::SparseMatrix<scalar_t, Eigen::RowMajor>::InnerIterator it(
304 : matrix, i);
305 656959 : it; ++it) {
306 621727 : eigenvectors.insert(it.row() + offset_rows, it.col() + offset_cols) =
307 621727 : it.value();
308 : }
309 : }
310 1315 : offset_rows += matrix.rows();
311 1315 : offset_cols += matrix.cols();
312 : }
313 619 : eigenvectors.makeCompressed();
314 :
315 619 : assert(
316 : eigenvectors.nonZeros() ==
317 : std::accumulate(non_zeros_per_inner_index.begin(), non_zeros_per_inner_index.end(), 0));
318 :
319 : // Get the combined eigenenergy matrix
320 619 : eigenenergies.reserve(Eigen::VectorXi::Constant(num_cols, 1));
321 619 : Eigen::Index offset = 0;
322 1934 : for (const auto &matrix : eigenenergies_blocks) {
323 1315 : set_task_status("Combining eigenenergies...");
324 31529 : for (int i = 0; i < matrix.size(); ++i) {
325 30214 : eigenenergies.insert(i + offset, i + offset) = matrix(i);
326 : }
327 1315 : offset += matrix.size();
328 : }
329 619 : eigenenergies.makeCompressed();
330 :
331 : // Fix phase ambiguity
332 619 : std::vector<scalar_t> map_col_to_max(num_cols, 0);
333 35851 : for (int row = 0; row < eigenvectors.outerSize(); ++row) {
334 35232 : set_task_status("Normalizing eigenvector phases...");
335 35232 : for (typename Eigen::SparseMatrix<scalar_t, Eigen::RowMajor>::InnerIterator it(
336 : eigenvectors, row);
337 656959 : it; ++it) {
338 621727 : if (std::abs(it.value()) > std::abs(map_col_to_max[it.col()])) {
339 141256 : map_col_to_max[it.col()] = it.value();
340 : }
341 : }
342 : }
343 :
344 619 : Eigen::SparseMatrix<scalar_t, Eigen::RowMajor> phase_matrix;
345 619 : phase_matrix.resize(num_cols, num_cols);
346 619 : phase_matrix.reserve(Eigen::VectorXi::Constant(num_cols, 1));
347 30833 : for (int i = 0; i < num_cols; ++i) {
348 30214 : phase_matrix.insert(i, i) = std::abs(map_col_to_max[i]) / map_col_to_max[i];
349 : }
350 619 : phase_matrix.makeCompressed();
351 :
352 619 : set_task_status("Applying eigenvector phases...");
353 619 : eigenvectors = eigenvectors * phase_matrix;
354 619 : }
355 :
356 : // Store the diagonalized hamiltonian
357 621 : matrix = eigenenergies;
358 621 : basis = basis->transformed(eigenvectors);
359 :
360 621 : hamiltonian_is_diagonal = true;
361 :
362 621 : return *this;
363 621 : }
364 :
365 : template <typename Derived>
366 2569 : bool System<Derived>::is_diagonal() const {
367 2569 : if (hamiltonian_requires_construction) {
368 39 : construct_hamiltonian();
369 39 : hamiltonian_requires_construction = false;
370 : }
371 :
372 2569 : if (!hamiltonian_is_diagonal) {
373 730 : real_t numerical_precision = 100 * std::numeric_limits<real_t>::epsilon();
374 :
375 16529 : for (int row = 0; row < matrix.outerSize(); ++row) {
376 16433 : for (typename Eigen::SparseMatrix<scalar_t, Eigen::RowMajor>::InnerIterator it(matrix,
377 : row);
378 32888 : it; ++it) {
379 17089 : if (it.row() != it.col() && std::abs(it.value()) > numerical_precision) {
380 634 : return false;
381 : }
382 : }
383 : }
384 :
385 96 : hamiltonian_is_diagonal = true;
386 : }
387 :
388 1935 : return true;
389 : }
390 :
391 : // Explicit instantiation
392 : template class System<SystemAtom<double>>;
393 : template class System<SystemAtom<std::complex<double>>>;
394 : template class System<SystemPair<double>>;
395 : template class System<SystemPair<std::complex<double>>>;
396 : } // namespace pairinteraction
|