1use scirs2_core::ndarray::{Array1, Array2, Array3, ArrayView1, ArrayView2, Axis};
8use scirs2_core::parallel_ops::*;
9use scirs2_core::Complex64;
10use serde::{Deserialize, Serialize};
11use std::collections::HashMap;
12
13use crate::error::{Result, SimulatorError};
14use crate::pauli::{PauliOperatorSum, PauliString};
15use crate::scirs2_integration::SciRS2Backend;
16use crate::scirs2_sparse::{
17 SciRS2SparseSolver, SparseEigenResult, SparseMatrix, SparseSolverConfig,
18};
19
20#[derive(Debug, Clone, Serialize, Deserialize)]
22pub struct SpectralAnalysisResult {
23 pub eigenvalues: Vec<f64>,
25 pub eigenvectors: Array2<Complex64>,
27 pub energy_gaps: Vec<f64>,
29 pub ground_state_energy: f64,
31 pub spectral_gap: f64,
33 pub participation_ratios: Vec<f64>,
35 pub entanglement_entropies: Vec<f64>,
37 pub spectral_stats: SpectralStatistics,
39}
40
41#[derive(Debug, Clone, Serialize, Deserialize)]
43pub struct PhaseTransitionResult {
44 pub parameters: Vec<f64>,
46 pub ground_state_energies: Vec<f64>,
48 pub energy_gaps: Vec<f64>,
50 pub order_parameters: Vec<f64>,
52 pub fidelity_susceptibility: Vec<f64>,
54 pub critical_points: Vec<f64>,
56 pub phase_boundaries: Vec<(f64, f64)>,
58}
59
60#[derive(Debug, Clone, Serialize, Deserialize)]
62pub struct SpectralDensityResult {
63 pub energy_grid: Vec<f64>,
65 pub density: Vec<f64>,
67 pub local_dos: Array2<f64>,
69 pub integrated_dos: Vec<f64>,
71 pub mobility_edges: Vec<f64>,
73}
74
75#[derive(Debug, Clone, Serialize, Deserialize)]
77pub struct EntanglementSpectrumResult {
78 pub eigenvalues: Vec<f64>,
80 pub entropy: f64,
82 pub renyi_entropies: HashMap<String, f64>,
84 pub entanglement_gap: f64,
86 pub participation_ratio: f64,
88 pub bipartition: Vec<usize>,
90}
91
92#[derive(Debug, Clone, Default, Serialize, Deserialize)]
94pub struct SpectralStatistics {
95 pub level_spacing_mean: f64,
97 pub level_spacing_std: f64,
98 pub nn_spacing_ratio: f64,
100 pub spectral_rigidity: f64,
102 pub number_variance: f64,
104 pub spectral_form_factor: Vec<Complex64>,
106 pub thouless_time: f64,
108}
109
110#[derive(Debug, Clone, Serialize, Deserialize)]
112pub struct BandStructureResult {
113 pub k_points: Vec<Array1<f64>>,
115 pub energy_bands: Array2<f64>,
117 pub band_gaps: Vec<f64>,
119 pub dos: SpectralDensityResult,
121 pub effective_masses: Vec<f64>,
123}
124
125#[derive(Debug, Clone)]
127pub struct SpectralConfig {
128 pub num_eigenvalues: usize,
130 pub which: String,
132 pub tolerance: f64,
134 pub max_iterations: usize,
136 pub parallel: bool,
138 pub energy_resolution: f64,
140 pub broadening: f64,
142}
143
144impl Default for SpectralConfig {
145 fn default() -> Self {
146 Self {
147 num_eigenvalues: 10,
148 which: "smallest".to_string(),
149 tolerance: 1e-10,
150 max_iterations: 1000,
151 parallel: true,
152 energy_resolution: 0.01,
153 broadening: 0.1,
154 }
155 }
156}
157
158pub struct SciRS2SpectralAnalyzer {
160 backend: Option<SciRS2Backend>,
162 config: SpectralConfig,
164 eigenvalue_cache: HashMap<String, SparseEigenResult>,
166}
167
168impl SciRS2SpectralAnalyzer {
169 pub fn new(config: SpectralConfig) -> Result<Self> {
171 Ok(Self {
172 backend: None,
173 config,
174 eigenvalue_cache: HashMap::new(),
175 })
176 }
177
178 pub fn with_backend(mut self) -> Result<Self> {
180 self.backend = Some(SciRS2Backend::new());
181 Ok(self)
182 }
183
184 pub fn analyze_spectrum(
186 &mut self,
187 hamiltonian: &SparseMatrix,
188 ) -> Result<SpectralAnalysisResult> {
189 if !hamiltonian.is_square() {
190 return Err(SimulatorError::InvalidInput(
191 "Hamiltonian must be square".to_string(),
192 ));
193 }
194
195 let eigen_result = self.compute_eigenvalues(hamiltonian)?;
197
198 let eigenvalues = eigen_result.eigenvalues;
199 let eigenvectors = eigen_result.eigenvectors;
200
201 let energy_gaps = self.calculate_energy_gaps(&eigenvalues);
203
204 let ground_state_energy = eigenvalues[0];
206 let spectral_gap = if eigenvalues.len() > 1 {
207 eigenvalues[1] - eigenvalues[0]
208 } else {
209 0.0
210 };
211
212 let participation_ratios = self.calculate_participation_ratios(&eigenvectors)?;
214
215 let entanglement_entropies = if hamiltonian.shape.0 <= 4096 {
217 self.calculate_entanglement_entropies(&eigenvectors)?
218 } else {
219 vec![0.0; eigenvalues.len()]
220 };
221
222 let spectral_stats = self.calculate_spectral_statistics(&eigenvalues)?;
224
225 Ok(SpectralAnalysisResult {
226 eigenvalues,
227 eigenvectors,
228 energy_gaps,
229 ground_state_energy,
230 spectral_gap,
231 participation_ratios,
232 entanglement_entropies,
233 spectral_stats,
234 })
235 }
236
237 pub fn analyze_phase_transition<F>(
239 &mut self,
240 hamiltonian_generator: F,
241 parameter_range: (f64, f64),
242 num_points: usize,
243 ) -> Result<PhaseTransitionResult>
244 where
245 F: Fn(f64) -> Result<SparseMatrix> + Sync + Send,
246 {
247 let parameters: Vec<f64> = (0..num_points)
248 .map(|i| {
249 parameter_range.0
250 + (parameter_range.1 - parameter_range.0) * i as f64 / (num_points - 1) as f64
251 })
252 .collect();
253
254 let results: Result<Vec<_>> = if self.config.parallel {
256 parameters
257 .par_iter()
258 .map(|¶m| {
259 let hamiltonian = hamiltonian_generator(param)?;
260 let mut analyzer = Self::new(self.config.clone())?;
261 let analysis = analyzer.analyze_spectrum(&hamiltonian)?;
262 Ok((
263 analysis.ground_state_energy,
264 analysis.spectral_gap,
265 analysis.eigenvectors.column(0).to_owned(),
266 ))
267 })
268 .collect()
269 } else {
270 parameters
271 .iter()
272 .map(|¶m| {
273 let hamiltonian = hamiltonian_generator(param)?;
274 let analysis = self.analyze_spectrum(&hamiltonian)?;
275 Ok((
276 analysis.ground_state_energy,
277 analysis.spectral_gap,
278 analysis.eigenvectors.column(0).to_owned(),
279 ))
280 })
281 .collect()
282 };
283
284 let results = results?;
285 let (ground_state_energies, energy_gaps, ground_states): (Vec<_>, Vec<_>, Vec<_>) =
286 results.into_iter().multiunzip();
287
288 let reference_state = &ground_states[0];
290 let order_parameters: Vec<f64> = ground_states
291 .iter()
292 .map(|state| {
293 let overlap: Complex64 = reference_state
294 .iter()
295 .zip(state.iter())
296 .map(|(&ref_amp, &state_amp)| ref_amp.conj() * state_amp)
297 .sum();
298 overlap.norm_sqr()
299 })
300 .collect();
301
302 let fidelity_susceptibility =
304 self.calculate_fidelity_susceptibility(&ground_states, ¶meters)?;
305
306 let critical_points =
308 self.detect_critical_points(&energy_gaps, &fidelity_susceptibility, ¶meters)?;
309
310 let phase_boundaries = self.identify_phase_boundaries(&order_parameters, ¶meters)?;
312
313 Ok(PhaseTransitionResult {
314 parameters,
315 ground_state_energies,
316 energy_gaps,
317 order_parameters,
318 fidelity_susceptibility,
319 critical_points,
320 phase_boundaries,
321 })
322 }
323
324 pub fn calculate_spectral_density(
326 &mut self,
327 hamiltonian: &SparseMatrix,
328 energy_range: (f64, f64),
329 ) -> Result<SpectralDensityResult> {
330 let num_points =
331 ((energy_range.1 - energy_range.0) / self.config.energy_resolution) as usize;
332 let energy_grid: Vec<f64> = (0..num_points)
333 .map(|i| {
334 energy_range.0
335 + (energy_range.1 - energy_range.0) * i as f64 / (num_points - 1) as f64
336 })
337 .collect();
338
339 let eigen_result = self.compute_eigenvalues(hamiltonian)?;
341 let eigenvalues = eigen_result.eigenvalues;
342
343 let density: Vec<f64> = energy_grid
345 .iter()
346 .map(|&energy| {
347 eigenvalues
348 .iter()
349 .map(|&eigenval| {
350 let diff = energy - eigenval;
351 (-diff * diff / (2.0 * self.config.broadening * self.config.broadening))
352 .exp()
353 })
354 .sum::<f64>()
355 / (self.config.broadening * (2.0 * std::f64::consts::PI).sqrt())
356 })
357 .collect();
358
359 let local_dos = Array2::zeros((hamiltonian.shape.0.min(100), num_points));
361
362 let mut integrated_dos = vec![0.0; num_points];
364 let mut cumulative = 0.0;
365 for (i, &d) in density.iter().enumerate() {
366 cumulative += d * self.config.energy_resolution;
367 integrated_dos[i] = cumulative;
368 }
369
370 let mobility_edges = Vec::new();
372
373 Ok(SpectralDensityResult {
374 energy_grid,
375 density,
376 local_dos,
377 integrated_dos,
378 mobility_edges,
379 })
380 }
381
382 pub fn calculate_entanglement_spectrum(
384 &mut self,
385 state: &Array1<Complex64>,
386 bipartition: &[usize],
387 ) -> Result<EntanglementSpectrumResult> {
388 let total_qubits = (state.len() as f64).log2() as usize;
389 if state.len() != 1 << total_qubits {
390 return Err(SimulatorError::InvalidInput(
391 "State vector length must be a power of 2".to_string(),
392 ));
393 }
394
395 let subsystem_size = bipartition.len();
396 let env_size = total_qubits - subsystem_size;
397
398 if subsystem_size == 0 || subsystem_size >= total_qubits {
399 return Err(SimulatorError::InvalidInput(
400 "Invalid bipartition specification".to_string(),
401 ));
402 }
403
404 let subsystem_dim = 1 << subsystem_size;
406 let env_dim = 1 << env_size;
407
408 let mut rho_reduced = Array2::zeros((subsystem_dim, subsystem_dim));
410
411 for i in 0..subsystem_dim {
412 for j in 0..subsystem_dim {
413 let mut element = Complex64::new(0.0, 0.0);
414 for k in 0..env_dim {
415 let full_i = self.encode_full_state(i, k, bipartition, total_qubits);
416 let full_j = self.encode_full_state(j, k, bipartition, total_qubits);
417 element += state[full_i].conj() * state[full_j];
418 }
419 rho_reduced[[i, j]] = element;
420 }
421 }
422
423 let eigenvalues = self.diagonalize_hermitian_matrix(&rho_reduced)?;
425
426 let mut nonzero_eigenvalues: Vec<f64> =
428 eigenvalues.into_iter().filter(|&x| x > 1e-15).collect();
429 nonzero_eigenvalues.sort_by(|a, b| b.partial_cmp(a).unwrap());
430
431 let entropy = -nonzero_eigenvalues.iter().map(|&p| p * p.ln()).sum::<f64>();
433
434 let mut renyi_entropies = HashMap::new();
436 for &n in &[2, 3, 4, 5] {
437 let renyi_n = if n == 1 {
438 entropy
439 } else {
440 (1.0 / (1.0 - n as f64))
441 * nonzero_eigenvalues
442 .iter()
443 .map(|&p| p.powf(n as f64))
444 .sum::<f64>()
445 .ln()
446 };
447 renyi_entropies.insert(format!("S_{n}"), renyi_n);
448 }
449
450 let entanglement_gap = if nonzero_eigenvalues.len() > 1 {
452 nonzero_eigenvalues[0] - nonzero_eigenvalues[1]
453 } else {
454 0.0
455 };
456
457 let participation_ratio = if nonzero_eigenvalues.is_empty() {
459 0.0
460 } else {
461 let sum_p2 = nonzero_eigenvalues.iter().map(|&p| p * p).sum::<f64>();
462 1.0 / sum_p2
463 };
464
465 Ok(EntanglementSpectrumResult {
466 eigenvalues: nonzero_eigenvalues,
467 entropy,
468 renyi_entropies,
469 entanglement_gap,
470 participation_ratio,
471 bipartition: bipartition.to_vec(),
472 })
473 }
474
475 pub fn calculate_band_structure<F>(
477 &mut self,
478 hamiltonian_generator: F,
479 k_path: &[Array1<f64>],
480 num_bands: usize,
481 ) -> Result<BandStructureResult>
482 where
483 F: Fn(&Array1<f64>) -> Result<SparseMatrix> + Sync + Send,
484 {
485 let band_results: Result<Vec<_>> = if self.config.parallel {
487 k_path
488 .par_iter()
489 .map(|k_point| {
490 let hamiltonian = hamiltonian_generator(k_point)?;
491 let mut config = SparseSolverConfig::default();
492 config.method = crate::scirs2_sparse::SparseSolverMethod::Lanczos;
493 let mut solver = SciRS2SparseSolver::new(config)?;
494 let result =
495 solver.solve_eigenvalue_problem(&hamiltonian, num_bands, "smallest")?;
496 Ok(result.eigenvalues)
497 })
498 .collect()
499 } else {
500 k_path
501 .iter()
502 .map(|k_point| {
503 let hamiltonian = hamiltonian_generator(k_point)?;
504 let mut config = SparseSolverConfig::default();
505 config.method = crate::scirs2_sparse::SparseSolverMethod::Lanczos;
506 let mut solver = SciRS2SparseSolver::new(config)?;
507 let result =
508 solver.solve_eigenvalue_problem(&hamiltonian, num_bands, "smallest")?;
509 Ok(result.eigenvalues)
510 })
511 .collect()
512 };
513
514 let band_results = band_results?;
515
516 let mut energy_bands = Array2::zeros((k_path.len(), num_bands));
518 for (i, eigenvalues) in band_results.iter().enumerate() {
519 for (j, &energy) in eigenvalues.iter().enumerate() {
520 if j < num_bands {
521 energy_bands[[i, j]] = energy;
522 }
523 }
524 }
525
526 let mut band_gaps = Vec::new();
528 for i in 0..num_bands - 1 {
529 let mut min_gap = f64::INFINITY;
530 for k in 0..k_path.len() {
531 let gap = energy_bands[[k, i + 1]] - energy_bands[[k, i]];
532 if gap < min_gap {
533 min_gap = gap;
534 }
535 }
536 band_gaps.push(min_gap);
537 }
538
539 let all_energies: Vec<f64> = energy_bands.iter().copied().collect();
541 let energy_range = (
542 all_energies.iter().copied().fold(f64::INFINITY, f64::min),
543 all_energies
544 .iter()
545 .copied()
546 .fold(f64::NEG_INFINITY, f64::max),
547 );
548
549 let num_energy_points = 200;
551 let energy_min = energy_range.0;
552 let energy_max = energy_range.1;
553 let energy_step = (energy_max - energy_min) / (num_energy_points - 1) as f64;
554
555 let mut energy_grid = Vec::new();
556 let mut density = Vec::new();
557
558 for i in 0..num_energy_points {
559 let energy = (i as f64).mul_add(energy_step, energy_min);
560 energy_grid.push(energy);
561
562 let window = energy_step * 2.0; let count = all_energies
565 .iter()
566 .filter(|&&e| (e - energy).abs() < window)
567 .count() as f64;
568
569 let dos_value = count / (window * k_path.len() as f64);
571 density.push(dos_value);
572 }
573
574 let mut integrated_dos = Vec::new();
576 let mut integral = 0.0;
577 for &dos_val in &density {
578 integral += dos_val * energy_step;
579 integrated_dos.push(integral);
580 }
581
582 let local_dos = Array2::zeros((all_energies.len().min(50), num_energy_points));
584
585 let dos = SpectralDensityResult {
586 energy_grid,
587 density,
588 local_dos,
589 integrated_dos,
590 mobility_edges: Vec::new(), };
592
593 let mut effective_masses = Vec::new();
595
596 for band in 0..num_bands {
597 let mut band_effective_masses = Vec::new();
598
599 for k_idx in 1..k_path.len() - 1 {
601 let e_prev = energy_bands[[k_idx - 1, band]];
602 let e_curr = energy_bands[[k_idx, band]];
603 let e_next = energy_bands[[k_idx + 1, band]];
604
605 let dk = 0.1; let second_derivative = (2.0f64.mul_add(-e_curr, e_next) + e_prev) / (dk * dk);
609
610 let effective_mass = if second_derivative.abs() > 1e-10 {
613 1.0 / second_derivative.abs()
614 } else {
615 f64::INFINITY };
617
618 band_effective_masses.push(effective_mass);
619 }
620
621 let avg_effective_mass = if band_effective_masses.is_empty() {
623 1.0 } else {
625 band_effective_masses.iter().sum::<f64>() / band_effective_masses.len() as f64
626 };
627
628 effective_masses.push(avg_effective_mass);
629 }
630
631 Ok(BandStructureResult {
632 k_points: k_path.to_vec(),
633 energy_bands,
634 band_gaps,
635 dos,
636 effective_masses,
637 })
638 }
639
640 fn compute_eigenvalues(&mut self, hamiltonian: &SparseMatrix) -> Result<SparseEigenResult> {
642 let cache_key = format!("eigenvals_{}_{}", hamiltonian.shape.0, hamiltonian.nnz);
643
644 if let Some(cached_result) = self.eigenvalue_cache.get(&cache_key) {
645 return Ok(cached_result.clone());
646 }
647
648 let mut config = SparseSolverConfig::default();
649 config.method = if hamiltonian.is_hermitian {
650 crate::scirs2_sparse::SparseSolverMethod::Lanczos
651 } else {
652 crate::scirs2_sparse::SparseSolverMethod::Arnoldi
653 };
654
655 let mut solver = if self.backend.is_some() {
656 SciRS2SparseSolver::new(config)?.with_backend()?
657 } else {
658 SciRS2SparseSolver::new(config)?
659 };
660
661 let result = solver.solve_eigenvalue_problem(
662 hamiltonian,
663 self.config.num_eigenvalues,
664 &self.config.which,
665 )?;
666
667 self.eigenvalue_cache.insert(cache_key, result.clone());
668 Ok(result)
669 }
670
671 fn calculate_energy_gaps(&self, eigenvalues: &[f64]) -> Vec<f64> {
673 eigenvalues
674 .windows(2)
675 .map(|window| window[1] - window[0])
676 .collect()
677 }
678
679 fn calculate_participation_ratios(&self, eigenvectors: &Array2<Complex64>) -> Result<Vec<f64>> {
681 let mut ratios = Vec::new();
682
683 for col_idx in 0..eigenvectors.ncols() {
684 let column = eigenvectors.column(col_idx);
685 let sum_p4: f64 = column.iter().map(|&c| c.norm_sqr().powi(2)).sum();
686 let ratio = if sum_p4 > 0.0 { 1.0 / sum_p4 } else { 0.0 };
687 ratios.push(ratio);
688 }
689
690 Ok(ratios)
691 }
692
693 fn calculate_entanglement_entropies(
695 &mut self,
696 eigenvectors: &Array2<Complex64>,
697 ) -> Result<Vec<f64>> {
698 let mut entropies = Vec::new();
699
700 for col_idx in 0..eigenvectors.ncols() {
701 let state = eigenvectors.column(col_idx).to_owned();
702 let num_qubits = (state.len() as f64).log2() as usize;
703
704 if num_qubits <= 10 {
705 let bipartition: Vec<usize> = (0..num_qubits / 2).collect();
707 let ent_result = self.calculate_entanglement_spectrum(&state, &bipartition);
708 let entropy = ent_result.map(|r| r.entropy).unwrap_or(0.0);
709 entropies.push(entropy);
710 } else {
711 entropies.push(0.0);
712 }
713 }
714
715 Ok(entropies)
716 }
717
718 fn calculate_spectral_statistics(&self, eigenvalues: &[f64]) -> Result<SpectralStatistics> {
720 if eigenvalues.len() < 3 {
721 return Ok(SpectralStatistics::default());
722 }
723
724 let spacings: Vec<f64> = eigenvalues
726 .windows(2)
727 .map(|window| window[1] - window[0])
728 .collect();
729
730 let level_spacing_mean = spacings.iter().sum::<f64>() / spacings.len() as f64;
731 let level_spacing_std = {
732 let variance = spacings
733 .iter()
734 .map(|&s| (s - level_spacing_mean).powi(2))
735 .sum::<f64>()
736 / spacings.len() as f64;
737 variance.sqrt()
738 };
739
740 let nn_spacing_ratio = if spacings.len() > 1 {
742 let ratios: Vec<f64> = spacings
743 .windows(2)
744 .map(|window| window[0].min(window[1]) / window[0].max(window[1]))
745 .collect();
746 ratios.iter().sum::<f64>() / ratios.len() as f64
747 } else {
748 0.0
749 };
750
751 let spectral_rigidity = level_spacing_std / level_spacing_mean;
753
754 let number_variance = eigenvalues.len() as f64 * level_spacing_std.powi(2);
756
757 let max_time_steps = 100;
759 let max_time = 10.0 * level_spacing_mean; let dt = max_time / max_time_steps as f64;
761
762 let mut spectral_form_factor = Vec::new();
763 let mut thouless_time = 1.0 / level_spacing_mean; for t_idx in 0..max_time_steps {
766 let t = t_idx as f64 * dt;
767
768 let mut sum_exp = Complex64::new(0.0, 0.0);
770 for &energy in eigenvalues {
771 let phase = energy * t;
772 sum_exp += Complex64::new(phase.cos(), phase.sin());
773 }
774
775 let form_factor_value = sum_exp.norm_sqr() / eigenvalues.len() as f64;
776 spectral_form_factor.push(Complex64::new(form_factor_value, 0.0));
777
778 if t_idx >= 5 && t_idx < spectral_form_factor.len().saturating_sub(5) {
780 let prev_values: Vec<f64> = spectral_form_factor[t_idx.saturating_sub(5)..t_idx]
782 .iter()
783 .map(|c| c.re)
784 .collect();
785 let next_values: Vec<f64> = spectral_form_factor
786 [t_idx..(t_idx + 5).min(spectral_form_factor.len())]
787 .iter()
788 .map(|c| c.re)
789 .collect();
790
791 let prev_slope =
793 (prev_values.last().unwrap() - prev_values.first().unwrap()) / (5.0 * dt);
794 let next_slope =
795 (next_values.last().unwrap() - next_values.first().unwrap()) / (5.0 * dt);
796
797 if prev_slope > 0.1 && next_slope.abs() < 0.05 && t > dt * 5.0 {
799 thouless_time = t;
800 break;
801 }
802 }
803 }
804
805 Ok(SpectralStatistics {
806 level_spacing_mean,
807 level_spacing_std,
808 nn_spacing_ratio,
809 spectral_rigidity,
810 number_variance,
811 spectral_form_factor,
812 thouless_time,
813 })
814 }
815
816 fn calculate_fidelity_susceptibility(
818 &self,
819 ground_states: &[Array1<Complex64>],
820 parameters: &[f64],
821 ) -> Result<Vec<f64>> {
822 let mut susceptibilities = Vec::new();
823
824 for i in 1..ground_states.len() - 1 {
825 let dparam = parameters[i + 1] - parameters[i - 1];
826 if dparam.abs() < 1e-15 {
827 susceptibilities.push(0.0);
828 continue;
829 }
830
831 let mut derivative = Array1::zeros(ground_states[i].len());
833 for j in 0..ground_states[i].len() {
834 derivative[j] = (ground_states[i + 1][j] - ground_states[i - 1][j]) / dparam;
835 }
836
837 let overlap: Complex64 = ground_states[i]
839 .iter()
840 .zip(derivative.iter())
841 .map(|(&gs, &deriv)| gs.conj() * deriv)
842 .sum();
843
844 let susceptibility =
845 derivative.iter().map(|&d| d.norm_sqr()).sum::<f64>() - overlap.norm_sqr();
846 susceptibilities.push(susceptibility);
847 }
848
849 if !susceptibilities.is_empty() {
851 susceptibilities.insert(0, susceptibilities[0]);
852 susceptibilities.push(susceptibilities[susceptibilities.len() - 1]);
853 }
854
855 Ok(susceptibilities)
856 }
857
858 fn detect_critical_points(
860 &self,
861 energy_gaps: &[f64],
862 fidelity_susceptibility: &[f64],
863 parameters: &[f64],
864 ) -> Result<Vec<f64>> {
865 let mut critical_points = Vec::new();
866
867 for i in 1..energy_gaps.len() - 1 {
869 if energy_gaps[i] < energy_gaps[i - 1] && energy_gaps[i] < energy_gaps[i + 1] {
870 critical_points.push(parameters[i]);
871 }
872 }
873
874 for i in 1..fidelity_susceptibility.len() - 1 {
876 if fidelity_susceptibility[i] > fidelity_susceptibility[i - 1]
877 && fidelity_susceptibility[i] > fidelity_susceptibility[i + 1]
878 {
879 critical_points.push(parameters[i]);
880 }
881 }
882
883 critical_points.sort_by(|a, b| a.partial_cmp(b).unwrap());
885 critical_points.dedup_by(|a, b| (*a - *b).abs() < 1e-6);
886
887 Ok(critical_points)
888 }
889
890 fn identify_phase_boundaries(
892 &self,
893 order_parameters: &[f64],
894 parameters: &[f64],
895 ) -> Result<Vec<(f64, f64)>> {
896 let mut boundaries = Vec::new();
897 let threshold = 0.1; for i in 1..order_parameters.len() {
900 let change = (order_parameters[i] - order_parameters[i - 1]).abs();
901 if change > threshold {
902 boundaries.push((parameters[i - 1], parameters[i]));
903 }
904 }
905
906 Ok(boundaries)
907 }
908
909 fn encode_full_state(
911 &self,
912 sub_idx: usize,
913 env_idx: usize,
914 bipartition: &[usize],
915 total_qubits: usize,
916 ) -> usize {
917 let mut full_idx = 0;
918 let mut sub_bit = 0;
919 let mut env_bit = 0;
920
921 for qubit in 0..total_qubits {
922 if bipartition.contains(&qubit) {
923 if (sub_idx >> sub_bit) & 1 == 1 {
924 full_idx |= 1 << (total_qubits - 1 - qubit);
925 }
926 sub_bit += 1;
927 } else {
928 if (env_idx >> env_bit) & 1 == 1 {
929 full_idx |= 1 << (total_qubits - 1 - qubit);
930 }
931 env_bit += 1;
932 }
933 }
934
935 full_idx
936 }
937
938 fn diagonalize_hermitian_matrix(&self, matrix: &Array2<Complex64>) -> Result<Vec<f64>> {
940 let n = matrix.nrows();
941 if n != matrix.ncols() {
942 return Err(SimulatorError::InvalidInput(
943 "Matrix must be square".to_string(),
944 ));
945 }
946
947 if n <= 8 {
949 let mut eigenvalues = Vec::new();
950
951 let mut x = Array1::from_vec(vec![Complex64::new(1.0, 0.0); n]);
953 let max_iterations = 100;
954 let tolerance = 1e-10;
955
956 for _ in 0..max_iterations {
957 let new_x = matrix.dot(&x);
959
960 let norm = new_x.iter().map(|c| c.norm_sqr()).sum::<f64>().sqrt();
962 if norm < tolerance {
963 break;
964 }
965
966 for (old, new) in x.iter_mut().zip(new_x.iter()) {
967 *old = *new / norm;
968 }
969 }
970
971 let ax = matrix.dot(&x);
973 let numerator: Complex64 = x
974 .iter()
975 .zip(ax.iter())
976 .map(|(xi, axi)| xi.conj() * axi)
977 .sum();
978 let denominator: f64 = x.iter().map(|xi| xi.norm_sqr()).sum();
979
980 if denominator > tolerance {
981 eigenvalues.push(numerator.re / denominator);
982 }
983
984 if n == 2 {
986 let trace = matrix[[0, 0]].re + matrix[[1, 1]].re;
987 let det = (matrix[[0, 0]] * matrix[[1, 1]] - matrix[[0, 1]] * matrix[[1, 0]]).re;
988
989 let discriminant = trace.mul_add(trace, -(4.0 * det));
991 if discriminant >= 0.0 {
992 let sqrt_disc = discriminant.sqrt();
993 eigenvalues.clear();
994 eigenvalues.push(f64::midpoint(trace, sqrt_disc));
995 eigenvalues.push((trace - sqrt_disc) / 2.0);
996 }
997 }
998
999 if eigenvalues.is_empty() {
1001 eprintln!(
1002 "Warning: Failed to compute eigenvalues properly, using diagonal approximation"
1003 );
1004 for i in 0..n {
1005 eigenvalues.push(matrix[[i, i]].re);
1006 }
1007 }
1008
1009 eigenvalues.sort_by(|a, b| b.partial_cmp(a).unwrap());
1010 Ok(eigenvalues)
1011 } else {
1012 Err(SimulatorError::UnsupportedOperation(
1014 format!("Matrix size {n} too large. Recommend using ndarray-linalg or LAPACK for proper eigenvalue computation")
1015 ))
1016 }
1017 }
1018
1019 pub const fn get_config(&self) -> &SpectralConfig {
1021 &self.config
1022 }
1023
1024 pub fn set_config(&mut self, config: SpectralConfig) {
1026 self.config = config;
1027 }
1028
1029 pub fn clear_cache(&mut self) {
1031 self.eigenvalue_cache.clear();
1032 }
1033}
1034
1035pub struct QuantumHamiltonianLibrary;
1037
1038impl QuantumHamiltonianLibrary {
1039 pub fn transverse_field_ising(num_sites: usize, j: f64, h: f64) -> Result<SparseMatrix> {
1041 use crate::scirs2_sparse::SparseMatrixUtils;
1042
1043 let mut pauli_terms = Vec::new();
1044
1045 for i in 0..num_sites - 1 {
1047 let mut pauli_string = "I".repeat(num_sites);
1048 pauli_string.replace_range(i..=i, "Z");
1049 pauli_string.replace_range(i + 1..i + 2, "Z");
1050 pauli_terms.push((pauli_string, -j));
1051 }
1052
1053 for i in 0..num_sites {
1055 let mut pauli_string = "I".repeat(num_sites);
1056 pauli_string.replace_range(i..=i, "X");
1057 pauli_terms.push((pauli_string, -h));
1058 }
1059
1060 SparseMatrixUtils::hamiltonian_from_pauli_strings(num_sites, &pauli_terms)
1061 }
1062
1063 pub fn heisenberg_model(num_sites: usize, jx: f64, jy: f64, jz: f64) -> Result<SparseMatrix> {
1065 use crate::scirs2_sparse::SparseMatrixUtils;
1066
1067 let mut pauli_terms = Vec::new();
1068
1069 for i in 0..num_sites - 1 {
1070 let mut xx_string = "I".repeat(num_sites);
1072 xx_string.replace_range(i..=i, "X");
1073 xx_string.replace_range(i + 1..i + 2, "X");
1074 pauli_terms.push((xx_string, jx));
1075
1076 let mut yy_string = "I".repeat(num_sites);
1078 yy_string.replace_range(i..=i, "Y");
1079 yy_string.replace_range(i + 1..i + 2, "Y");
1080 pauli_terms.push((yy_string, jy));
1081
1082 let mut zz_string = "I".repeat(num_sites);
1084 zz_string.replace_range(i..=i, "Z");
1085 zz_string.replace_range(i + 1..i + 2, "Z");
1086 pauli_terms.push((zz_string, jz));
1087 }
1088
1089 SparseMatrixUtils::hamiltonian_from_pauli_strings(num_sites, &pauli_terms)
1090 }
1091
1092 pub fn random_matrix_model(size: usize, density: f64, hermitian: bool) -> SparseMatrix {
1094 use crate::scirs2_sparse::SparseMatrixUtils;
1095 SparseMatrixUtils::random_sparse(size, density, hermitian)
1096 }
1097}
1098
1099trait MultiUnzip<A, B, C> {
1101 fn multiunzip(self) -> (Vec<A>, Vec<B>, Vec<C>);
1102}
1103
1104impl<I, A, B, C> MultiUnzip<A, B, C> for I
1105where
1106 I: Iterator<Item = (A, B, C)>,
1107{
1108 fn multiunzip(self) -> (Vec<A>, Vec<B>, Vec<C>) {
1109 let mut vec_a = Vec::new();
1110 let mut vec_b = Vec::new();
1111 let mut vec_c = Vec::new();
1112
1113 for (a, b, c) in self {
1114 vec_a.push(a);
1115 vec_b.push(b);
1116 vec_c.push(c);
1117 }
1118
1119 (vec_a, vec_b, vec_c)
1120 }
1121}
1122
1123pub fn benchmark_spectral_analysis(system_size: usize) -> Result<HashMap<String, f64>> {
1125 let mut results = HashMap::new();
1126
1127 let hamiltonian = QuantumHamiltonianLibrary::transverse_field_ising(system_size, 1.0, 0.5)?;
1129
1130 let config = SpectralConfig {
1132 num_eigenvalues: 10.min(hamiltonian.shape.0),
1133 ..Default::default()
1134 };
1135
1136 let mut analyzer = SciRS2SpectralAnalyzer::new(config)?;
1137
1138 let start_time = std::time::Instant::now();
1140 let _analysis = analyzer.analyze_spectrum(&hamiltonian)?;
1141 let analysis_time = start_time.elapsed().as_secs_f64();
1142
1143 results.insert("SpectrumAnalysis".to_string(), analysis_time);
1144
1145 let start_time = std::time::Instant::now();
1147 let _density = analyzer.calculate_spectral_density(&hamiltonian, (-5.0, 5.0))?;
1148 let density_time = start_time.elapsed().as_secs_f64();
1149
1150 results.insert("SpectralDensity".to_string(), density_time);
1151
1152 Ok(results)
1153}
1154
1155#[cfg(test)]
1156mod tests {
1157 use super::*;
1158 use approx::assert_abs_diff_eq;
1159
1160 #[test]
1161 fn test_spectral_analyzer_creation() {
1162 let config = SpectralConfig::default();
1163 let analyzer = SciRS2SpectralAnalyzer::new(config).unwrap();
1164 assert_eq!(analyzer.config.num_eigenvalues, 10);
1165 }
1166
1167 #[test]
1168 fn test_transverse_field_ising() {
1169 let hamiltonian = QuantumHamiltonianLibrary::transverse_field_ising(3, 1.0, 0.5).unwrap();
1170 assert_eq!(hamiltonian.shape, (8, 8));
1171 assert!(hamiltonian.is_hermitian);
1172 }
1173
1174 #[test]
1175 fn test_heisenberg_model() {
1176 let hamiltonian = QuantumHamiltonianLibrary::heisenberg_model(2, 1.0, 1.0, 1.0).unwrap();
1177 assert_eq!(hamiltonian.shape, (4, 4));
1178 assert!(hamiltonian.is_hermitian);
1179 }
1180
1181 #[test]
1182 fn test_spectral_analysis() {
1183 let hamiltonian = QuantumHamiltonianLibrary::transverse_field_ising(3, 1.0, 0.5).unwrap();
1184 let config = SpectralConfig {
1185 num_eigenvalues: 5,
1186 ..Default::default()
1187 };
1188
1189 let mut analyzer = SciRS2SpectralAnalyzer::new(config).unwrap();
1190 let result = analyzer.analyze_spectrum(&hamiltonian).unwrap();
1191
1192 assert_eq!(result.eigenvalues.len(), 5);
1193 assert!(result.spectral_gap >= 0.0);
1194 assert_eq!(result.participation_ratios.len(), 5);
1195 }
1196
1197 #[test]
1198 fn test_entanglement_spectrum() {
1199 let mut state = Array1::zeros(4);
1200 state[0] = Complex64::new(1.0 / 2.0_f64.sqrt(), 0.0);
1201 state[3] = Complex64::new(1.0 / 2.0_f64.sqrt(), 0.0);
1202
1203 let config = SpectralConfig::default();
1204 let mut analyzer = SciRS2SpectralAnalyzer::new(config).unwrap();
1205
1206 let bipartition = vec![0];
1207 let result = analyzer
1208 .calculate_entanglement_spectrum(&state, &bipartition)
1209 .unwrap();
1210
1211 assert!(result.entropy > 0.0);
1212 assert_eq!(result.bipartition, vec![0]);
1213 }
1214
1215 #[test]
1216 fn test_energy_gaps() {
1217 let eigenvalues = vec![0.0, 1.0, 3.0, 6.0];
1218 let config = SpectralConfig::default();
1219 let analyzer = SciRS2SpectralAnalyzer::new(config).unwrap();
1220
1221 let gaps = analyzer.calculate_energy_gaps(&eigenvalues);
1222
1223 assert_eq!(gaps, vec![1.0, 2.0, 3.0]);
1224 }
1225
1226 #[test]
1227 fn test_participation_ratios() {
1228 let mut eigenvectors = Array2::zeros((4, 2));
1229 eigenvectors[[0, 0]] = Complex64::new(1.0, 0.0);
1230 eigenvectors[[1, 1]] = Complex64::new(0.5, 0.0);
1231 eigenvectors[[2, 1]] = Complex64::new(0.5, 0.0);
1232 eigenvectors[[3, 1]] = Complex64::new(0.5, 0.0);
1233
1234 let config = SpectralConfig::default();
1235 let analyzer = SciRS2SpectralAnalyzer::new(config).unwrap();
1236
1237 let ratios = analyzer
1238 .calculate_participation_ratios(&eigenvectors)
1239 .unwrap();
1240
1241 assert_eq!(ratios.len(), 2);
1242 assert_abs_diff_eq!(ratios[0], 1.0, epsilon = 1e-10); assert!(ratios[1] > 1.0); }
1245}