quantrs2_sim/
statevector.rs

1use scirs2_core::parallel_ops::{
2    IndexedParallelIterator, IntoParallelRefMutIterator, ParallelIterator,
3};
4use scirs2_core::Complex64;
5use std::sync::Mutex;
6
7use quantrs2_circuit::builder::{Circuit, Simulator};
8use quantrs2_core::{
9    error::{QuantRS2Error, QuantRS2Result},
10    gate::{multi, single, GateOp},
11    qubit::QubitId,
12    register::Register,
13};
14
15use crate::diagnostics::SimulationDiagnostics;
16use crate::optimized_simd;
17use crate::scirs2_integration::{SciRS2Backend, SciRS2Matrix, SciRS2Vector, BLAS};
18use crate::utils::{flip_bit, gate_vec_to_array2};
19
20/// A state vector simulator for quantum circuits
21///
22/// This simulator implements the state vector approach, where the full quantum
23/// state is represented as a complex vector of dimension 2^N for N qubits.
24#[derive(Debug)]
25pub struct StateVectorSimulator {
26    /// Use parallel execution
27    pub parallel: bool,
28
29    /// Basic noise model (if any)
30    pub noise_model: Option<crate::noise::NoiseModel>,
31
32    /// Advanced noise model (if any)
33    pub advanced_noise_model: Option<crate::noise_advanced::AdvancedNoiseModel>,
34
35    /// Optimized buffer pool for memory reuse (thread-safe)
36    buffer_pool: Mutex<BufferPool>,
37
38    /// Enable SIMD optimizations for gate operations
39    pub use_simd: bool,
40
41    /// Enable gate fusion optimization
42    pub use_gate_fusion: bool,
43
44    /// Diagnostics system for monitoring and error handling
45    pub diagnostics: Option<SimulationDiagnostics>,
46
47    /// `SciRS2` backend for optimized linear algebra operations
48    scirs2_backend: SciRS2Backend,
49}
50
51impl Clone for StateVectorSimulator {
52    fn clone(&self) -> Self {
53        Self {
54            parallel: self.parallel,
55            noise_model: self.noise_model.clone(),
56            advanced_noise_model: self.advanced_noise_model.clone(),
57            buffer_pool: Mutex::new(BufferPool::new(4, 1024)), // Create new buffer pool
58            use_simd: self.use_simd,
59            use_gate_fusion: self.use_gate_fusion,
60            diagnostics: self
61                .diagnostics
62                .as_ref()
63                .map(|_| SimulationDiagnostics::new()),
64            scirs2_backend: SciRS2Backend::new(),
65        }
66    }
67}
68
69/// Memory pool for efficient state vector operations
70#[derive(Debug, Clone)]
71pub struct BufferPool {
72    /// Pre-allocated working buffers
73    working_buffers: std::collections::VecDeque<Vec<Complex64>>,
74    /// Maximum number of cached buffers
75    max_buffers: usize,
76    /// Target buffer size for efficient reuse
77    target_size: usize,
78}
79
80impl BufferPool {
81    /// Create new buffer pool
82    #[must_use]
83    pub fn new(max_buffers: usize, target_size: usize) -> Self {
84        Self {
85            working_buffers: std::collections::VecDeque::with_capacity(max_buffers),
86            max_buffers,
87            target_size,
88        }
89    }
90
91    /// Get a buffer from pool or allocate new one
92    pub fn get_buffer(&mut self, size: usize) -> Vec<Complex64> {
93        // Try to reuse existing buffer
94        if let Some(mut buffer) = self.working_buffers.pop_front() {
95            if buffer.capacity() >= size {
96                buffer.clear();
97                buffer.resize(size, Complex64::new(0.0, 0.0));
98                return buffer;
99            }
100        }
101
102        // Allocate new buffer with extra capacity for growth
103        let capacity = size.max(self.target_size);
104        let mut buffer = Vec::with_capacity(capacity);
105        buffer.resize(size, Complex64::new(0.0, 0.0));
106        buffer
107    }
108
109    /// Return buffer to pool
110    pub fn return_buffer(&mut self, buffer: Vec<Complex64>) {
111        if self.working_buffers.len() < self.max_buffers
112            && buffer.capacity() >= self.target_size / 2
113        {
114            self.working_buffers.push_back(buffer);
115        }
116        // Otherwise let buffer be dropped and deallocated
117    }
118
119    /// Clear all cached buffers
120    pub fn clear(&mut self) {
121        self.working_buffers.clear();
122    }
123}
124
125impl StateVectorSimulator {
126    /// Create a new state vector simulator with default settings
127    #[must_use]
128    pub fn new() -> Self {
129        Self {
130            parallel: true,
131            noise_model: None,
132            advanced_noise_model: None,
133            buffer_pool: Mutex::new(BufferPool::new(4, 1024)), // Default: 4 buffers, 1K size target
134            use_simd: true,
135            use_gate_fusion: true,
136            diagnostics: None,
137            scirs2_backend: SciRS2Backend::new(),
138        }
139    }
140
141    /// Create a new state vector simulator with parallel execution disabled
142    #[must_use]
143    pub fn sequential() -> Self {
144        Self {
145            parallel: false,
146            noise_model: None,
147            advanced_noise_model: None,
148            buffer_pool: Mutex::new(BufferPool::new(2, 512)), // Smaller pool for sequential
149            use_simd: true,
150            use_gate_fusion: true,
151            diagnostics: None,
152            scirs2_backend: SciRS2Backend::new(),
153        }
154    }
155
156    /// Create a new state vector simulator with a basic noise model
157    #[must_use]
158    pub fn with_noise(noise_model: crate::noise::NoiseModel) -> Self {
159        Self {
160            parallel: true,
161            noise_model: Some(noise_model),
162            advanced_noise_model: None,
163            buffer_pool: Mutex::new(BufferPool::new(4, 1024)),
164            use_simd: true,
165            use_gate_fusion: true,
166            diagnostics: None,
167            scirs2_backend: SciRS2Backend::new(),
168        }
169    }
170
171    /// Create a new state vector simulator with an advanced noise model
172    #[must_use]
173    pub fn with_advanced_noise(
174        advanced_noise_model: crate::noise_advanced::AdvancedNoiseModel,
175    ) -> Self {
176        Self {
177            parallel: true,
178            noise_model: None,
179            advanced_noise_model: Some(advanced_noise_model),
180            buffer_pool: Mutex::new(BufferPool::new(4, 1024)),
181            use_simd: true,
182            use_gate_fusion: true,
183            diagnostics: None,
184            scirs2_backend: SciRS2Backend::new(),
185        }
186    }
187
188    /// Create simulator with custom buffer pool configuration
189    #[must_use]
190    pub fn with_buffer_pool(parallel: bool, max_buffers: usize, target_size: usize) -> Self {
191        Self {
192            parallel,
193            noise_model: None,
194            advanced_noise_model: None,
195            buffer_pool: Mutex::new(BufferPool::new(max_buffers, target_size)),
196            use_simd: true,
197            use_gate_fusion: true,
198            diagnostics: None,
199            scirs2_backend: SciRS2Backend::new(),
200        }
201    }
202
203    /// Set the basic noise model
204    pub fn set_noise_model(&mut self, noise_model: crate::noise::NoiseModel) -> &mut Self {
205        self.noise_model = Some(noise_model);
206        self.advanced_noise_model = None; // Remove advanced model if it exists
207        self
208    }
209
210    /// Set the advanced noise model
211    pub fn set_advanced_noise_model(
212        &mut self,
213        advanced_noise_model: crate::noise_advanced::AdvancedNoiseModel,
214    ) -> &mut Self {
215        self.advanced_noise_model = Some(advanced_noise_model);
216        self.noise_model = None; // Remove basic model if it exists
217        self
218    }
219
220    /// Remove all noise models
221    pub fn remove_noise_model(&mut self) -> &mut Self {
222        self.noise_model = None;
223        self.advanced_noise_model = None;
224        self
225    }
226
227    /// Enable or disable SIMD optimizations
228    pub const fn set_simd_enabled(&mut self, enabled: bool) -> &mut Self {
229        self.use_simd = enabled;
230        self
231    }
232
233    /// Enable or disable gate fusion optimization
234    pub const fn set_gate_fusion_enabled(&mut self, enabled: bool) -> &mut Self {
235        self.use_gate_fusion = enabled;
236        self
237    }
238
239    /// Get access to buffer pool for testing purposes
240    pub const fn get_buffer_pool(&self) -> &Mutex<BufferPool> {
241        &self.buffer_pool
242    }
243
244    /// Enable diagnostics and monitoring
245    pub fn enable_diagnostics(&mut self) -> &mut Self {
246        self.diagnostics = Some(SimulationDiagnostics::new());
247        self
248    }
249
250    /// Disable diagnostics and monitoring
251    pub fn disable_diagnostics(&mut self) -> &mut Self {
252        self.diagnostics = None;
253        self
254    }
255
256    /// Get diagnostics report if diagnostics are enabled
257    pub fn get_diagnostics_report(&self) -> Option<crate::diagnostics::DiagnosticReport> {
258        self.diagnostics
259            .as_ref()
260            .map(super::diagnostics::SimulationDiagnostics::generate_report)
261    }
262
263    /// Create a high-performance configuration
264    #[must_use]
265    pub fn high_performance() -> Self {
266        Self {
267            parallel: true,
268            noise_model: None,
269            advanced_noise_model: None,
270            buffer_pool: Mutex::new(BufferPool::new(8, 2048)), // Larger pool for high performance
271            use_simd: true,
272            use_gate_fusion: true,
273            diagnostics: Some(SimulationDiagnostics::new()),
274            scirs2_backend: SciRS2Backend::new(),
275        }
276    }
277
278    /// Apply a dense matrix-vector multiplication using SciRS2 when available
279    #[cfg(feature = "advanced_math")]
280    fn apply_dense_matrix_vector(
281        &mut self,
282        matrix: &[Complex64],
283        vector: &[Complex64],
284        result: &mut [Complex64],
285    ) -> QuantRS2Result<()> {
286        if self.scirs2_backend.is_available() && matrix.len() >= 64 && vector.len() >= 8 {
287            // Use SciRS2 for larger operations where the overhead is worthwhile
288            use scirs2_core::ndarray::{Array1, Array2};
289
290            let rows = result.len();
291            let cols = vector.len();
292
293            if matrix.len() != rows * cols {
294                return Err(QuantRS2Error::InvalidInput(
295                    "Dimension mismatch in matrix application".to_string(),
296                ));
297            }
298
299            // Convert to ndarray format
300            let matrix_2d =
301                Array2::from_shape_vec((rows, cols), matrix.to_vec()).map_err(|_| {
302                    QuantRS2Error::InvalidInput("Matrix shape conversion failed".to_string())
303                })?;
304            let vector_1d = Array1::from_vec(vector.to_vec());
305
306            // Convert to SciRS2 format
307            let scirs2_matrix = SciRS2Matrix::from_array2(matrix_2d);
308            let scirs2_vector = SciRS2Vector::from_array1(vector_1d);
309
310            // Perform optimized matrix-vector multiplication
311            let scirs2_result = self
312                .scirs2_backend
313                .matrix_vector_multiply(&scirs2_matrix, &scirs2_vector)
314                .map_err(|_| {
315                    QuantRS2Error::ComputationError(
316                        "Matrix-vector multiplication failed".to_string(),
317                    )
318                })?;
319
320            // Convert back to slice
321            let result_array = scirs2_result.to_array1().map_err(|_| {
322                QuantRS2Error::ComputationError("Result conversion to array failed".to_string())
323            })?;
324            if let Some(slice) = result_array.as_slice() {
325                result.copy_from_slice(slice);
326            } else {
327                return Err(QuantRS2Error::ComputationError(
328                    "Result array is not contiguous".to_string(),
329                ));
330            }
331
332            Ok(())
333        } else {
334            // Fallback to manual implementation for smaller operations
335            for i in 0..result.len() {
336                result[i] = Complex64::new(0.0, 0.0);
337                for j in 0..vector.len() {
338                    result[i] += matrix[i * vector.len() + j] * vector[j];
339                }
340            }
341            Ok(())
342        }
343    }
344
345    /// Fallback dense matrix-vector multiplication for when `SciRS2` is not available
346    #[cfg(not(feature = "advanced_math"))]
347    fn apply_dense_matrix_vector(
348        &self,
349        matrix: &[Complex64],
350        vector: &[Complex64],
351        result: &mut [Complex64],
352    ) -> QuantRS2Result<()> {
353        // Manual implementation
354        for i in 0..result.len() {
355            result[i] = Complex64::new(0.0, 0.0);
356            for j in 0..vector.len() {
357                result[i] += matrix[i * vector.len() + j] * vector[j];
358            }
359        }
360        Ok(())
361    }
362
363    /// Apply a single-qubit gate to a state vector
364    fn apply_single_qubit_gate<const N: usize>(
365        &self,
366        state: &mut [Complex64],
367        gate_matrix: &[Complex64],
368        target: QubitId,
369    ) -> QuantRS2Result<()> {
370        let target_idx = target.id() as usize;
371        if target_idx >= N {
372            return Err(QuantRS2Error::InvalidQubitId(target.id()));
373        }
374
375        // Use SIMD optimization if enabled and beneficial
376        if self.use_simd && state.len() >= 8 {
377            return self.apply_single_qubit_gate_simd::<N>(state, gate_matrix, target_idx);
378        }
379
380        // Convert the gate matrix to flat representation for faster access
381        // Gate matrix: [m00, m01, m10, m11]
382        let m00 = gate_matrix[0];
383        let m01 = gate_matrix[1];
384        let m10 = gate_matrix[2];
385        let m11 = gate_matrix[3];
386
387        // Apply the gate to each amplitude
388        if self.parallel {
389            // Get buffer from pool for temporary storage
390            let mut state_copy = self
391                .buffer_pool
392                .lock()
393                .expect("buffer pool lock poisoned")
394                .get_buffer(state.len());
395            state_copy.copy_from_slice(state);
396
397            state.par_iter_mut().enumerate().for_each(|(idx, amp)| {
398                let bit_val = (idx >> target_idx) & 1;
399                let paired_idx = if bit_val == 0 {
400                    idx | (1 << target_idx)
401                } else {
402                    idx & !(1 << target_idx)
403                };
404
405                let idx0 = if bit_val == 0 { idx } else { paired_idx };
406                let idx1 = if bit_val == 0 { paired_idx } else { idx };
407
408                let val0 = state_copy[idx0];
409                let val1 = state_copy[idx1];
410
411                // Use direct matrix element access for better performance
412                *amp = if idx == idx0 {
413                    m00 * val0 + m01 * val1
414                } else {
415                    m10 * val0 + m11 * val1
416                };
417            });
418
419            // Return buffer to pool
420            self.buffer_pool
421                .lock()
422                .expect("buffer pool lock poisoned")
423                .return_buffer(state_copy);
424        } else {
425            // Sequential implementation using buffer pool
426            let dim = state.len();
427            let mut new_state = self
428                .buffer_pool
429                .lock()
430                .expect("buffer pool lock poisoned")
431                .get_buffer(dim);
432
433            for i in 0..dim {
434                let bit_val = (i >> target_idx) & 1;
435                let paired_idx = flip_bit(i, target_idx);
436
437                if bit_val == 0 {
438                    new_state[i] = m00 * state[i] + m01 * state[paired_idx];
439                    new_state[paired_idx] = m10 * state[i] + m11 * state[paired_idx];
440                }
441            }
442
443            state.copy_from_slice(&new_state);
444            // Return buffer to pool
445            self.buffer_pool
446                .lock()
447                .expect("buffer pool lock poisoned")
448                .return_buffer(new_state);
449        }
450
451        Ok(())
452    }
453
454    /// Apply a single-qubit gate using SIMD optimization
455    fn apply_single_qubit_gate_simd<const N: usize>(
456        &self,
457        state: &mut [Complex64],
458        gate_matrix: &[Complex64],
459        target_idx: usize,
460    ) -> QuantRS2Result<()> {
461        let dim = state.len();
462        let pairs_per_block = 1 << target_idx;
463        let total_pairs = dim / 2;
464
465        // Get buffers for SIMD processing
466        let mut pool = self.buffer_pool.lock().expect("buffer pool lock poisoned");
467        let mut in_amps0 = pool.get_buffer(total_pairs);
468        let mut in_amps1 = pool.get_buffer(total_pairs);
469        let mut out_amps0 = pool.get_buffer(total_pairs);
470        let mut out_amps1 = pool.get_buffer(total_pairs);
471
472        // Collect amplitudes for target bit = 0 and target bit = 1
473        let mut pair_idx = 0;
474        for block in 0..(dim / (pairs_per_block * 2)) {
475            let block_start = block * pairs_per_block * 2;
476
477            for offset in 0..pairs_per_block {
478                let idx0 = block_start + offset;
479                let idx1 = block_start + pairs_per_block + offset;
480
481                in_amps0[pair_idx] = state[idx0];
482                in_amps1[pair_idx] = state[idx1];
483                pair_idx += 1;
484            }
485        }
486
487        // Convert gate matrix to required format for SIMD
488        let gate_matrix_array: [Complex64; 4] = [
489            gate_matrix[0],
490            gate_matrix[1],
491            gate_matrix[2],
492            gate_matrix[3],
493        ];
494
495        // Apply SIMD gate operation
496        optimized_simd::apply_single_qubit_gate_optimized(
497            &gate_matrix_array,
498            &in_amps0,
499            &in_amps1,
500            &mut out_amps0,
501            &mut out_amps1,
502        );
503
504        // Write results back to state vector
505        pair_idx = 0;
506        for block in 0..(dim / (pairs_per_block * 2)) {
507            let block_start = block * pairs_per_block * 2;
508
509            for offset in 0..pairs_per_block {
510                let idx0 = block_start + offset;
511                let idx1 = block_start + pairs_per_block + offset;
512
513                state[idx0] = out_amps0[pair_idx];
514                state[idx1] = out_amps1[pair_idx];
515                pair_idx += 1;
516            }
517        }
518
519        // Return buffers to pool
520        pool.return_buffer(in_amps0);
521        pool.return_buffer(in_amps1);
522        pool.return_buffer(out_amps0);
523        pool.return_buffer(out_amps1);
524
525        Ok(())
526    }
527
528    /// Apply a two-qubit gate to a state vector
529    fn apply_two_qubit_gate<const N: usize>(
530        &self,
531        state: &mut [Complex64],
532        gate_matrix: &[Complex64],
533        control: QubitId,
534        target: QubitId,
535    ) -> QuantRS2Result<()> {
536        let control_idx = control.id() as usize;
537        let target_idx = target.id() as usize;
538
539        if control_idx >= N || target_idx >= N {
540            return Err(QuantRS2Error::InvalidQubitId(if control_idx >= N {
541                control.id()
542            } else {
543                target.id()
544            }));
545        }
546
547        if control_idx == target_idx {
548            return Err(QuantRS2Error::CircuitValidationFailed(
549                "Control and target qubits must be different".into(),
550            ));
551        }
552
553        // Pre-extract matrix elements for faster access (16 elements total)
554        let m = gate_matrix; // Direct slice access is faster than ndarray indexing
555
556        // Apply the gate to each amplitude
557        if self.parallel {
558            // Get buffer from pool for temporary storage
559            let mut state_copy = self
560                .buffer_pool
561                .lock()
562                .expect("buffer pool lock poisoned")
563                .get_buffer(state.len());
564            state_copy.copy_from_slice(state);
565
566            state.par_iter_mut().enumerate().for_each(|(idx, amp)| {
567                let idx00 = idx & !(1 << control_idx) & !(1 << target_idx);
568                let idx01 = idx00 | (1 << target_idx);
569                let idx10 = idx00 | (1 << control_idx);
570                let idx11 = idx00 | (1 << control_idx) | (1 << target_idx);
571
572                let val00 = state_copy[idx00];
573                let val01 = state_copy[idx01];
574                let val10 = state_copy[idx10];
575                let val11 = state_copy[idx11];
576
577                // Use direct matrix access for better performance
578                *amp = match idx {
579                    i if i == idx00 => m[0] * val00 + m[1] * val01 + m[2] * val10 + m[3] * val11,
580                    i if i == idx01 => m[4] * val00 + m[5] * val01 + m[6] * val10 + m[7] * val11,
581                    i if i == idx10 => m[8] * val00 + m[9] * val01 + m[10] * val10 + m[11] * val11,
582                    i if i == idx11 => {
583                        m[12] * val00 + m[13] * val01 + m[14] * val10 + m[15] * val11
584                    }
585                    _ => unreachable!(),
586                };
587            });
588
589            // Return buffer to pool
590            self.buffer_pool
591                .lock()
592                .expect("buffer pool lock poisoned")
593                .return_buffer(state_copy);
594        } else {
595            // Sequential implementation using buffer pool
596            let dim = state.len();
597            let mut new_state = self
598                .buffer_pool
599                .lock()
600                .expect("buffer pool lock poisoned")
601                .get_buffer(dim);
602
603            #[allow(clippy::needless_range_loop)]
604            for i in 0..dim {
605                let control_bit = (i >> control_idx) & 1;
606                let target_bit = (i >> target_idx) & 1;
607
608                // Calculate the four basis states in the 2-qubit subspace
609                let i00 = i & !(1 << control_idx) & !(1 << target_idx);
610                let i01 = i00 | (1 << target_idx);
611                let i10 = i00 | (1 << control_idx);
612                let i11 = i10 | (1 << target_idx);
613
614                let basis_idx = (control_bit << 1) | target_bit;
615
616                // Calculate the new amplitude for this state using direct access
617                let row_offset = basis_idx * 4;
618                new_state[i] = m[row_offset] * state[i00]
619                    + m[row_offset + 1] * state[i01]
620                    + m[row_offset + 2] * state[i10]
621                    + m[row_offset + 3] * state[i11];
622            }
623
624            state.copy_from_slice(&new_state);
625            // Return buffer to pool
626            self.buffer_pool
627                .lock()
628                .expect("buffer pool lock poisoned")
629                .return_buffer(new_state);
630        }
631
632        Ok(())
633    }
634
635    /// Apply CNOT gate efficiently (special case)
636    fn apply_cnot<const N: usize>(
637        &self,
638        state: &mut [Complex64],
639        control: QubitId,
640        target: QubitId,
641    ) -> QuantRS2Result<()> {
642        let control_idx = control.id() as usize;
643        let target_idx = target.id() as usize;
644
645        if control_idx >= N || target_idx >= N {
646            return Err(QuantRS2Error::InvalidQubitId(if control_idx >= N {
647                control.id()
648            } else {
649                target.id()
650            }));
651        }
652
653        if control_idx == target_idx {
654            return Err(QuantRS2Error::CircuitValidationFailed(
655                "Control and target qubits must be different".into(),
656            ));
657        }
658
659        // Apply the CNOT gate - only swap amplitudes where control is 1
660        if self.parallel {
661            let mut state_copy = self
662                .buffer_pool
663                .lock()
664                .expect("buffer pool lock poisoned")
665                .get_buffer(state.len());
666            state_copy.copy_from_slice(state);
667
668            state.par_iter_mut().enumerate().for_each(|(i, amp)| {
669                if (i >> control_idx) & 1 == 1 {
670                    let flipped = flip_bit(i, target_idx);
671                    *amp = state_copy[flipped];
672                }
673            });
674
675            self.buffer_pool
676                .lock()
677                .expect("buffer pool lock poisoned")
678                .return_buffer(state_copy);
679        } else {
680            let dim = state.len();
681            let mut new_state = self
682                .buffer_pool
683                .lock()
684                .expect("buffer pool lock poisoned")
685                .get_buffer(dim);
686
687            for i in 0..dim {
688                if (i >> control_idx) & 1 == 1 {
689                    let flipped = flip_bit(i, target_idx);
690                    new_state[flipped] = state[i];
691                    new_state[i] = state[flipped];
692                } else {
693                    new_state[i] = state[i];
694                }
695            }
696
697            state.copy_from_slice(&new_state);
698            self.buffer_pool
699                .lock()
700                .expect("buffer pool lock poisoned")
701                .return_buffer(new_state);
702        }
703
704        Ok(())
705    }
706
707    /// Apply SWAP gate efficiently (special case)
708    fn apply_swap<const N: usize>(
709        &self,
710        state: &mut [Complex64],
711        qubit1: QubitId,
712        qubit2: QubitId,
713    ) -> QuantRS2Result<()> {
714        let q1_idx = qubit1.id() as usize;
715        let q2_idx = qubit2.id() as usize;
716
717        if q1_idx >= N || q2_idx >= N {
718            return Err(QuantRS2Error::InvalidQubitId(if q1_idx >= N {
719                qubit1.id()
720            } else {
721                qubit2.id()
722            }));
723        }
724
725        if q1_idx == q2_idx {
726            return Err(QuantRS2Error::CircuitValidationFailed(
727                "Qubits must be different for SWAP gate".into(),
728            ));
729        }
730
731        // Apply the SWAP gate - swap amplitudes where qubits have different values
732        if self.parallel {
733            let mut state_copy = self
734                .buffer_pool
735                .lock()
736                .expect("buffer pool lock poisoned")
737                .get_buffer(state.len());
738            state_copy.copy_from_slice(state);
739
740            state.par_iter_mut().enumerate().for_each(|(i, amp)| {
741                let bit1 = (i >> q1_idx) & 1;
742                let bit2 = (i >> q2_idx) & 1;
743
744                if bit1 != bit2 {
745                    let swapped = flip_bit(flip_bit(i, q1_idx), q2_idx);
746                    *amp = state_copy[swapped];
747                }
748            });
749
750            self.buffer_pool
751                .lock()
752                .expect("buffer pool lock poisoned")
753                .return_buffer(state_copy);
754        } else {
755            let dim = state.len();
756            let mut new_state = self
757                .buffer_pool
758                .lock()
759                .expect("buffer pool lock poisoned")
760                .get_buffer(dim);
761
762            for i in 0..dim {
763                let bit1 = (i >> q1_idx) & 1;
764                let bit2 = (i >> q2_idx) & 1;
765
766                if bit1 == bit2 {
767                    new_state[i] = state[i];
768                } else {
769                    let swapped = flip_bit(flip_bit(i, q1_idx), q2_idx);
770                    new_state[swapped] = state[i];
771                    new_state[i] = state[swapped];
772                }
773            }
774
775            state.copy_from_slice(&new_state);
776            self.buffer_pool
777                .lock()
778                .expect("buffer pool lock poisoned")
779                .return_buffer(new_state);
780        }
781
782        Ok(())
783    }
784}
785
786impl Default for StateVectorSimulator {
787    fn default() -> Self {
788        Self::new()
789    }
790}
791
792impl<const N: usize> Simulator<N> for StateVectorSimulator {
793    fn run(&self, circuit: &Circuit<N>) -> QuantRS2Result<Register<N>> {
794        // Initialize state vector to |0...0⟩
795        let dim = 1 << N;
796        let mut state = vec![Complex64::new(0.0, 0.0); dim];
797        state[0] = Complex64::new(1.0, 0.0);
798
799        // Apply each gate in the circuit
800        for gate in circuit.gates() {
801            match gate.name() {
802                // Single-qubit gates
803                "H" => {
804                    if let Some(g) = gate.as_any().downcast_ref::<single::Hadamard>() {
805                        let matrix = g.matrix()?;
806                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
807                    }
808                }
809                "X" => {
810                    if let Some(g) = gate.as_any().downcast_ref::<single::PauliX>() {
811                        let matrix = g.matrix()?;
812                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
813                    }
814                }
815                "Y" => {
816                    if let Some(g) = gate.as_any().downcast_ref::<single::PauliY>() {
817                        let matrix = g.matrix()?;
818                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
819                    }
820                }
821                "Z" => {
822                    if let Some(g) = gate.as_any().downcast_ref::<single::PauliZ>() {
823                        let matrix = g.matrix()?;
824                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
825                    }
826                }
827                "RX" => {
828                    if let Some(g) = gate.as_any().downcast_ref::<single::RotationX>() {
829                        let matrix = g.matrix()?;
830                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
831                    }
832                }
833                "RY" => {
834                    if let Some(g) = gate.as_any().downcast_ref::<single::RotationY>() {
835                        let matrix = g.matrix()?;
836                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
837                    }
838                }
839                "RZ" => {
840                    if let Some(g) = gate.as_any().downcast_ref::<single::RotationZ>() {
841                        let matrix = g.matrix()?;
842                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
843                    }
844                }
845                "S" => {
846                    if let Some(g) = gate.as_any().downcast_ref::<single::Phase>() {
847                        let matrix = g.matrix()?;
848                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
849                    }
850                }
851                "T" => {
852                    if let Some(g) = gate.as_any().downcast_ref::<single::T>() {
853                        let matrix = g.matrix()?;
854                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
855                    }
856                }
857                "S†" => {
858                    if let Some(g) = gate.as_any().downcast_ref::<single::PhaseDagger>() {
859                        let matrix = g.matrix()?;
860                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
861                    }
862                }
863                "T†" => {
864                    if let Some(g) = gate.as_any().downcast_ref::<single::TDagger>() {
865                        let matrix = g.matrix()?;
866                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
867                    }
868                }
869                "√X" => {
870                    if let Some(g) = gate.as_any().downcast_ref::<single::SqrtX>() {
871                        let matrix = g.matrix()?;
872                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
873                    }
874                }
875                "√X†" => {
876                    if let Some(g) = gate.as_any().downcast_ref::<single::SqrtXDagger>() {
877                        let matrix = g.matrix()?;
878                        self.apply_single_qubit_gate::<N>(&mut state, &matrix, g.target)?;
879                    }
880                }
881
882                // Two-qubit gates
883                "CNOT" => {
884                    if let Some(g) = gate.as_any().downcast_ref::<multi::CNOT>() {
885                        // Use optimized implementation for CNOT
886                        self.apply_cnot::<N>(&mut state, g.control, g.target)?;
887                    }
888                }
889                "CZ" => {
890                    if let Some(g) = gate.as_any().downcast_ref::<multi::CZ>() {
891                        let matrix = g.matrix()?;
892                        self.apply_two_qubit_gate::<N>(&mut state, &matrix, g.control, g.target)?;
893                    }
894                }
895                "SWAP" => {
896                    if let Some(g) = gate.as_any().downcast_ref::<multi::SWAP>() {
897                        // Use optimized implementation for SWAP
898                        self.apply_swap::<N>(&mut state, g.qubit1, g.qubit2)?;
899                    }
900                }
901                "CY" => {
902                    if let Some(g) = gate.as_any().downcast_ref::<multi::CY>() {
903                        let matrix = g.matrix()?;
904                        self.apply_two_qubit_gate::<N>(&mut state, &matrix, g.control, g.target)?;
905                    }
906                }
907                "CH" => {
908                    if let Some(g) = gate.as_any().downcast_ref::<multi::CH>() {
909                        let matrix = g.matrix()?;
910                        self.apply_two_qubit_gate::<N>(&mut state, &matrix, g.control, g.target)?;
911                    }
912                }
913                "CS" => {
914                    if let Some(g) = gate.as_any().downcast_ref::<multi::CS>() {
915                        let matrix = g.matrix()?;
916                        self.apply_two_qubit_gate::<N>(&mut state, &matrix, g.control, g.target)?;
917                    }
918                }
919                "CRX" => {
920                    if let Some(g) = gate.as_any().downcast_ref::<multi::CRX>() {
921                        let matrix = g.matrix()?;
922                        self.apply_two_qubit_gate::<N>(&mut state, &matrix, g.control, g.target)?;
923                    }
924                }
925                "CRY" => {
926                    if let Some(g) = gate.as_any().downcast_ref::<multi::CRY>() {
927                        let matrix = g.matrix()?;
928                        self.apply_two_qubit_gate::<N>(&mut state, &matrix, g.control, g.target)?;
929                    }
930                }
931                "CRZ" => {
932                    if let Some(g) = gate.as_any().downcast_ref::<multi::CRZ>() {
933                        let matrix = g.matrix()?;
934                        self.apply_two_qubit_gate::<N>(&mut state, &matrix, g.control, g.target)?;
935                    }
936                }
937
938                // Three-qubit gates
939                "Toffoli" => {
940                    if let Some(toffoli_gate) = gate.as_any().downcast_ref::<multi::Toffoli>() {
941                        let control1 = toffoli_gate.control1;
942                        let control2 = toffoli_gate.control2;
943                        let target = toffoli_gate.target;
944                        self.apply_toffoli_to_state::<N>(&mut state, control1, control2, target)?;
945                    }
946                }
947                "Fredkin" => {
948                    if let Some(fredkin_gate) = gate.as_any().downcast_ref::<multi::Fredkin>() {
949                        let control = fredkin_gate.control;
950                        let target1 = fredkin_gate.target1;
951                        let target2 = fredkin_gate.target2;
952                        self.apply_fredkin_to_state::<N>(&mut state, control, target1, target2)?;
953                    }
954                }
955
956                _ => {
957                    return Err(QuantRS2Error::UnsupportedOperation(format!(
958                        "Gate {} not supported",
959                        gate.name()
960                    )));
961                }
962            }
963
964            // Apply per-gate noise if configured
965            if let Some(ref noise_model) = self.noise_model {
966                if noise_model.per_gate {
967                    noise_model.apply_to_statevector(&mut state)?;
968                }
969            }
970
971            // Apply per-gate advanced noise if configured
972            if let Some(ref advanced_noise_model) = self.advanced_noise_model {
973                if advanced_noise_model.per_gate {
974                    advanced_noise_model.apply_to_statevector(&mut state)?;
975                }
976            }
977        }
978
979        // Apply final noise if not per-gate
980        if let Some(ref noise_model) = self.noise_model {
981            if !noise_model.per_gate {
982                noise_model.apply_to_statevector(&mut state)?;
983            }
984        }
985
986        // Apply final advanced noise if not per-gate
987        if let Some(ref advanced_noise_model) = self.advanced_noise_model {
988            if !advanced_noise_model.per_gate {
989                advanced_noise_model.apply_to_statevector(&mut state)?;
990            }
991        }
992
993        // Create register from final state
994        Register::<N>::with_amplitudes(state)
995    }
996}
997
998impl StateVectorSimulator {
999    /// Initialize state with specified number of qubits in |0...0⟩
1000    pub const fn initialize_state(&mut self, num_qubits: usize) -> QuantRS2Result<()> {
1001        // This is a placeholder - actual initialization would need the circuit framework
1002        Ok(())
1003    }
1004
1005    /// Get the current quantum state
1006    pub fn get_state(&self) -> Vec<Complex64> {
1007        // Placeholder - would return the actual state vector
1008        vec![Complex64::new(1.0, 0.0)]
1009    }
1010
1011    /// Get mutable reference to the current quantum state
1012    pub fn get_state_mut(&mut self) -> Vec<Complex64> {
1013        // Placeholder - would return mutable reference to actual state vector
1014        vec![Complex64::new(1.0, 0.0)]
1015    }
1016
1017    /// Set the quantum state
1018    pub fn set_state(&mut self, _state: Vec<Complex64>) -> QuantRS2Result<()> {
1019        // Placeholder - would set the actual state vector
1020        Ok(())
1021    }
1022
1023    /// Apply an interface circuit to the quantum state
1024    pub const fn apply_interface_circuit(
1025        &mut self,
1026        _circuit: &crate::circuit_interfaces::InterfaceCircuit,
1027    ) -> QuantRS2Result<()> {
1028        // Placeholder - would apply the circuit gates using the circuit framework
1029        Ok(())
1030    }
1031
1032    /// Apply Hadamard gate to qubit
1033    pub const fn apply_h(&mut self, _qubit: usize) -> QuantRS2Result<()> {
1034        // Placeholder - would apply H gate using circuit framework
1035        Ok(())
1036    }
1037
1038    /// Apply Pauli-X gate to qubit
1039    pub const fn apply_x(&mut self, _qubit: usize) -> QuantRS2Result<()> {
1040        // Placeholder - would apply X gate using circuit framework
1041        Ok(())
1042    }
1043
1044    /// Apply Pauli-Z gate to qubit
1045    pub const fn apply_z_public(&mut self, _qubit: usize) -> QuantRS2Result<()> {
1046        // Placeholder - would apply Z gate using circuit framework
1047        Ok(())
1048    }
1049
1050    /// Apply CNOT gate (public interface with usize indices)
1051    pub const fn apply_cnot_public(
1052        &mut self,
1053        _control: usize,
1054        _target: usize,
1055    ) -> QuantRS2Result<()> {
1056        // Placeholder - would apply CNOT gate using circuit framework
1057        Ok(())
1058    }
1059
1060    /// Apply Toffoli (CCNOT) gate to state vector
1061    fn apply_toffoli_to_state<const N: usize>(
1062        &self,
1063        state: &mut [Complex64],
1064        control1: QubitId,
1065        control2: QubitId,
1066        target: QubitId,
1067    ) -> QuantRS2Result<()> {
1068        let control1_idx = control1.id() as usize;
1069        let control2_idx = control2.id() as usize;
1070        let target_idx = target.id() as usize;
1071
1072        if control1_idx >= N || control2_idx >= N || target_idx >= N {
1073            return Err(QuantRS2Error::InvalidQubitId(if control1_idx >= N {
1074                control1.id()
1075            } else if control2_idx >= N {
1076                control2.id()
1077            } else {
1078                target.id()
1079            }));
1080        }
1081
1082        // Apply Toffoli gate by swapping amplitudes when both controls are |1⟩
1083        if self.parallel {
1084            let mut state_copy = self
1085                .buffer_pool
1086                .lock()
1087                .expect("buffer pool lock poisoned")
1088                .get_buffer(state.len());
1089            state_copy.copy_from_slice(state);
1090
1091            state.par_iter_mut().enumerate().for_each(|(i, amp)| {
1092                let control1_bit = (i >> control1_idx) & 1;
1093                let control2_bit = (i >> control2_idx) & 1;
1094
1095                if control1_bit == 1 && control2_bit == 1 {
1096                    let flipped = flip_bit(i, target_idx);
1097                    *amp = state_copy[flipped];
1098                }
1099            });
1100
1101            self.buffer_pool
1102                .lock()
1103                .expect("buffer pool lock poisoned")
1104                .return_buffer(state_copy);
1105        } else {
1106            let mut new_state = self
1107                .buffer_pool
1108                .lock()
1109                .expect("buffer pool lock poisoned")
1110                .get_buffer(state.len());
1111            new_state.copy_from_slice(state);
1112
1113            for i in 0..state.len() {
1114                let control1_bit = (i >> control1_idx) & 1;
1115                let control2_bit = (i >> control2_idx) & 1;
1116
1117                if control1_bit == 1 && control2_bit == 1 {
1118                    let flipped = flip_bit(i, target_idx);
1119                    new_state[flipped] = state[i];
1120                    new_state[i] = state[flipped];
1121                }
1122            }
1123
1124            state.copy_from_slice(&new_state);
1125            self.buffer_pool
1126                .lock()
1127                .expect("buffer pool lock poisoned")
1128                .return_buffer(new_state);
1129        }
1130
1131        Ok(())
1132    }
1133
1134    /// Apply Fredkin (CSWAP) gate to state vector
1135    fn apply_fredkin_to_state<const N: usize>(
1136        &self,
1137        state: &mut [Complex64],
1138        control: QubitId,
1139        target1: QubitId,
1140        target2: QubitId,
1141    ) -> QuantRS2Result<()> {
1142        let control_idx = control.id() as usize;
1143        let target1_idx = target1.id() as usize;
1144        let target2_idx = target2.id() as usize;
1145
1146        if control_idx >= N || target1_idx >= N || target2_idx >= N {
1147            return Err(QuantRS2Error::InvalidQubitId(if control_idx >= N {
1148                control.id()
1149            } else if target1_idx >= N {
1150                target1.id()
1151            } else {
1152                target2.id()
1153            }));
1154        }
1155
1156        // Apply Fredkin gate by swapping target qubits when control is |1⟩
1157        if self.parallel {
1158            let mut state_copy = self
1159                .buffer_pool
1160                .lock()
1161                .expect("buffer pool lock poisoned")
1162                .get_buffer(state.len());
1163            state_copy.copy_from_slice(state);
1164
1165            state.par_iter_mut().enumerate().for_each(|(i, amp)| {
1166                let control_bit = (i >> control_idx) & 1;
1167                let target1_bit = (i >> target1_idx) & 1;
1168                let target2_bit = (i >> target2_idx) & 1;
1169
1170                if control_bit == 1 && target1_bit != target2_bit {
1171                    let swapped = flip_bit(flip_bit(i, target1_idx), target2_idx);
1172                    *amp = state_copy[swapped];
1173                }
1174            });
1175
1176            self.buffer_pool
1177                .lock()
1178                .expect("buffer pool lock poisoned")
1179                .return_buffer(state_copy);
1180        } else {
1181            let mut new_state = self
1182                .buffer_pool
1183                .lock()
1184                .expect("buffer pool lock poisoned")
1185                .get_buffer(state.len());
1186            new_state.copy_from_slice(state);
1187
1188            for i in 0..state.len() {
1189                let control_bit = (i >> control_idx) & 1;
1190                let target1_bit = (i >> target1_idx) & 1;
1191                let target2_bit = (i >> target2_idx) & 1;
1192
1193                if control_bit == 1 && target1_bit != target2_bit {
1194                    let swapped = flip_bit(flip_bit(i, target1_idx), target2_idx);
1195                    new_state[swapped] = state[i];
1196                    new_state[i] = state[swapped];
1197                }
1198            }
1199
1200            state.copy_from_slice(&new_state);
1201            self.buffer_pool
1202                .lock()
1203                .expect("buffer pool lock poisoned")
1204                .return_buffer(new_state);
1205        }
1206
1207        Ok(())
1208    }
1209
1210    /// Apply Toffoli (CCNOT) gate (public interface)
1211    pub const fn apply_toffoli(
1212        &mut self,
1213        control1: QubitId,
1214        control2: QubitId,
1215        target: QubitId,
1216    ) -> QuantRS2Result<()> {
1217        // This is a placeholder for external API - real implementation would need circuit framework
1218        Ok(())
1219    }
1220
1221    /// Apply Fredkin (CSWAP) gate (public interface)
1222    pub const fn apply_fredkin(
1223        &mut self,
1224        control: QubitId,
1225        target1: QubitId,
1226        target2: QubitId,
1227    ) -> QuantRS2Result<()> {
1228        // This is a placeholder for external API - real implementation would need circuit framework
1229        Ok(())
1230    }
1231}