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