scirs2-optimize 0.4.2

Optimization module for SciRS2 (scirs2-optimize)
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
//! Clustering algorithms for local minima in multi-start optimization
//!
//! This module provides algorithms to identify, cluster, and analyze local minima
//! found during multi-start optimization strategies. It helps distinguish between
//! unique local minima and provides insights into the optimization landscape.

use crate::error::OptimizeError;
use crate::unconstrained::{minimize, Method, Options};
use scirs2_core::ndarray::{Array1, Array2, ArrayView1};
use std::collections::HashMap;

/// Configuration for clustering local minima
#[derive(Debug, Clone)]
pub struct ClusteringOptions {
    /// Distance threshold for considering minima as belonging to the same cluster
    pub distance_threshold: f64,
    /// Relative tolerance for function values when clustering
    pub function_tolerance: f64,
    /// Maximum number of clusters to form
    pub max_clusters: Option<usize>,
    /// Clustering algorithm to use
    pub algorithm: ClusteringAlgorithm,
    /// Whether to normalize coordinates before clustering
    pub normalize_coordinates: bool,
    /// Whether to use function values in clustering
    pub use_function_values: bool,
    /// Weight for function values vs coordinates in distance calculation
    pub function_weight: f64,
}

impl Default for ClusteringOptions {
    fn default() -> Self {
        Self {
            distance_threshold: 1e-3,
            function_tolerance: 1e-6,
            max_clusters: None,
            algorithm: ClusteringAlgorithm::Hierarchical,
            normalize_coordinates: true,
            use_function_values: true,
            function_weight: 0.1,
        }
    }
}

/// Clustering algorithms available
#[derive(Debug, Clone, Copy)]
pub enum ClusteringAlgorithm {
    /// Hierarchical clustering with single linkage
    Hierarchical,
    /// K-means clustering
    KMeans,
    /// Density-based clustering (DBSCAN-like)
    Density,
    /// Custom threshold-based clustering
    Threshold,
}

/// Represents a local minimum found during optimization
#[derive(Debug, Clone)]
pub struct LocalMinimum<S> {
    /// Location of the minimum
    pub x: Array1<f64>,
    /// Function value at the minimum
    pub f: f64,
    /// Original function value (for generic type)
    pub fun_value: S,
    /// Number of optimization iterations to reach this minimum
    pub nit: usize,
    /// Number of function evaluations
    pub func_evals: usize,
    /// Whether optimization was successful
    pub success: bool,
    /// Starting point that led to this minimum
    pub start_point: Array1<f64>,
    /// Cluster ID (assigned after clustering)
    pub cluster_id: Option<usize>,
    /// Distance to cluster centroid
    pub cluster_distance: Option<f64>,
}

/// Result of clustering analysis
#[derive(Debug, Clone)]
pub struct ClusteringResult<S> {
    /// All local minima found
    pub minima: Vec<LocalMinimum<S>>,
    /// Cluster centroids
    pub centroids: Vec<ClusterCentroid>,
    /// Number of clusters formed
    pub num_clusters: usize,
    /// Silhouette score (quality measure)
    pub silhouette_score: Option<f64>,
    /// Within-cluster sum of squares
    pub wcss: f64,
    /// Best minimum found (global optimum candidate)
    pub best_minimum: Option<LocalMinimum<S>>,
}

/// Cluster centroid information
#[derive(Debug, Clone)]
pub struct ClusterCentroid {
    /// Centroid coordinates
    pub x: Array1<f64>,
    /// Average function value in cluster
    pub f_avg: f64,
    /// Best (minimum) function value in cluster
    pub f_min: f64,
    /// Number of minima in this cluster
    pub size: usize,
    /// Cluster radius (max distance from centroid)
    pub radius: f64,
}

/// Multi-start optimization with clustering
#[allow(dead_code)]
pub fn multi_start_with_clustering<F, S>(
    fun: F,
    start_points: &[Array1<f64>],
    method: Method,
    options: Option<Options>,
    clustering_options: Option<ClusteringOptions>,
) -> Result<ClusteringResult<S>, OptimizeError>
where
    F: FnMut(&ArrayView1<f64>) -> S + Clone,
    S: Into<f64> + Clone + From<f64>,
{
    let clustering_opts = clustering_options.unwrap_or_default();
    let mut minima = Vec::new();

    // Run optimization from each starting point
    for start_point in start_points {
        let fun_clone = fun.clone();

        match minimize(
            fun_clone,
            start_point.as_slice().expect("Operation failed"),
            method,
            options.clone(),
        ) {
            Ok(result) => {
                let minimum = LocalMinimum {
                    x: result.x.clone(),
                    f: result.fun.clone().into(),
                    fun_value: result.fun,
                    nit: result.nit,
                    func_evals: result.func_evals,
                    success: result.success,
                    start_point: start_point.clone(),
                    cluster_id: None,
                    cluster_distance: None,
                };
                minima.push(minimum);
            }
            Err(_) => {
                // Skip failed optimizations but could log them
                continue;
            }
        }
    }

    if minima.is_empty() {
        return Err(OptimizeError::ComputationError(
            "No successful optimizations found".to_string(),
        ));
    }

    // Perform clustering
    cluster_minima(&mut minima, &clustering_opts)?;

    // Compute cluster centroids and statistics
    let centroids = compute_cluster_centroids(&minima)?;
    let num_clusters = centroids.len();
    let wcss = compute_wcss(&minima, &centroids);
    let silhouette_score = compute_silhouette_score(&minima);

    // Find best minimum
    let best_minimum = minima
        .iter()
        .filter(|m| m.success)
        .min_by(|a, b| a.f.partial_cmp(&b.f).expect("Operation failed"))
        .cloned();

    Ok(ClusteringResult {
        minima,
        centroids,
        num_clusters,
        silhouette_score,
        wcss,
        best_minimum,
    })
}

/// Cluster local minima using the specified algorithm
#[allow(dead_code)]
fn cluster_minima<S>(
    minima: &mut [LocalMinimum<S>],
    options: &ClusteringOptions,
) -> Result<(), OptimizeError>
where
    S: Clone,
{
    if minima.is_empty() {
        return Ok(());
    }

    match options.algorithm {
        ClusteringAlgorithm::Hierarchical => hierarchical_clustering(minima, options),
        ClusteringAlgorithm::KMeans => kmeans_clustering(minima, options),
        ClusteringAlgorithm::Density => density_clustering(minima, options),
        ClusteringAlgorithm::Threshold => threshold_clustering(minima, options),
    }
}

/// Hierarchical clustering implementation
#[allow(dead_code)]
fn hierarchical_clustering<S>(
    minima: &mut [LocalMinimum<S>],
    options: &ClusteringOptions,
) -> Result<(), OptimizeError>
where
    S: Clone,
{
    let n = minima.len();
    if n <= 1 {
        if n == 1 {
            minima[0].cluster_id = Some(0);
        }
        return Ok(());
    }

    // Compute distance matrix
    let distances = compute_distance_matrix(minima, options);

    // Perform hierarchical clustering using single linkage
    let mut cluster_assignments = vec![None; n];
    let mut next_cluster_id = n;

    // Initialize: each point is its own cluster
    for (i, assignment) in cluster_assignments.iter_mut().enumerate().take(n) {
        *assignment = Some(i);
    }

    // Build hierarchy by merging closest clusters
    let mut active_clusters: Vec<Vec<usize>> = (0..n).map(|i| vec![i]).collect();

    while active_clusters.len() > 1 {
        let mut min_dist = f64::INFINITY;
        let mut merge_pair = (0, 1);

        // Find closest pair of clusters
        for i in 0..active_clusters.len() {
            for j in (i + 1)..active_clusters.len() {
                let cluster_dist =
                    compute_cluster_distance(&active_clusters[i], &active_clusters[j], &distances);
                if cluster_dist < min_dist {
                    min_dist = cluster_dist;
                    merge_pair = (i, j);
                }
            }
        }

        // Stop if distance threshold exceeded
        if min_dist > options.distance_threshold {
            break;
        }

        // Check max clusters constraint
        if let Some(max_clusters) = options.max_clusters {
            if active_clusters.len() <= max_clusters {
                break;
            }
        }

        // Merge clusters
        let (i, j) = merge_pair;
        let mut merged_cluster = active_clusters[i].clone();
        merged_cluster.extend(&active_clusters[j]);

        // Update cluster assignments
        for &point_idx in &merged_cluster {
            cluster_assignments[point_idx] = Some(next_cluster_id);
        }

        // Remove old clusters and add merged one
        let mut new_clusters = Vec::new();
        for (k, cluster) in active_clusters.iter().enumerate() {
            if k != i && k != j {
                new_clusters.push(cluster.clone());
            }
        }
        new_clusters.push(merged_cluster);
        active_clusters = new_clusters;
        next_cluster_id += 1;
    }

    // Assign final cluster IDs (renumber to be sequential)
    let mut cluster_map = HashMap::new();
    let mut final_cluster_id = 0;

    for cluster_id in cluster_assignments.iter().flatten() {
        if !cluster_map.contains_key(cluster_id) {
            cluster_map.insert(*cluster_id, final_cluster_id);
            final_cluster_id += 1;
        }
    }

    // Update minima with final cluster assignments
    for (i, minimum) in minima.iter_mut().enumerate() {
        if let Some(cluster_id) = cluster_assignments[i] {
            minimum.cluster_id = cluster_map.get(&cluster_id).copied();
        }
    }

    Ok(())
}

/// K-means clustering implementation
#[allow(dead_code)]
fn kmeans_clustering<S>(
    minima: &mut [LocalMinimum<S>],
    options: &ClusteringOptions,
) -> Result<(), OptimizeError>
where
    S: Clone,
{
    let n = minima.len();
    if n <= 1 {
        if n == 1 {
            minima[0].cluster_id = Some(0);
        }
        return Ok(());
    }

    // Determine number of clusters
    let k = if let Some(max_k) = options.max_clusters {
        std::cmp::min(max_k, n)
    } else {
        // Use elbow method or default to sqrt(n)
        std::cmp::min((n as f64).sqrt() as usize + 1, n)
    };

    if k >= n {
        // Each point is its own cluster
        for (i, minimum) in minima.iter_mut().enumerate() {
            minimum.cluster_id = Some(i);
        }
        return Ok(());
    }

    // Get feature vectors for clustering
    let features = extract_features(minima, options);
    let dim = features.ncols();

    // Initialize centroids randomly (or use k-means++)
    let mut centroids = initialize_centroids_plus_plus(&features, k);
    let mut assignments = vec![0; n];
    let max_iter = 100;
    let tolerance = 1e-6;

    for _iter in 0..max_iter {
        let mut changed = false;

        // Assign points to nearest centroids
        for (i, assignment) in assignments.iter_mut().enumerate().take(n) {
            let mut min_dist = f64::INFINITY;
            let mut best_cluster = 0;

            for j in 0..k {
                let dist = euclidean_distance(&features.row(i), &centroids.row(j));
                if dist < min_dist {
                    min_dist = dist;
                    best_cluster = j;
                }
            }

            if *assignment != best_cluster {
                *assignment = best_cluster;
                changed = true;
            }
        }

        if !changed {
            break;
        }

        // Update centroids
        let mut new_centroids = Array2::zeros((k, dim));
        let mut cluster_sizes = vec![0; k];

        for i in 0..n {
            let cluster = assignments[i];
            cluster_sizes[cluster] += 1;
            for d in 0..dim {
                new_centroids[[cluster, d]] += features[[i, d]];
            }
        }

        for j in 0..k {
            if cluster_sizes[j] > 0 {
                for d in 0..dim {
                    new_centroids[[j, d]] /= cluster_sizes[j] as f64;
                }
            }
        }

        // Check convergence
        let centroid_change = (&centroids - &new_centroids).mapv(|x| x.abs()).sum();

        centroids = new_centroids;

        if centroid_change < tolerance {
            break;
        }
    }

    // Assign cluster IDs to minima
    for (i, minimum) in minima.iter_mut().enumerate() {
        minimum.cluster_id = Some(assignments[i]);
    }

    Ok(())
}

/// Density-based clustering (DBSCAN-like)
#[allow(dead_code)]
fn density_clustering<S>(
    minima: &mut [LocalMinimum<S>],
    options: &ClusteringOptions,
) -> Result<(), OptimizeError>
where
    S: Clone,
{
    let n = minima.len();
    if n <= 1 {
        if n == 1 {
            minima[0].cluster_id = Some(0);
        }
        return Ok(());
    }

    let eps = options.distance_threshold;
    let min_pts = 2; // Minimum points to form a cluster

    let distances = compute_distance_matrix(minima, options);
    let mut cluster_assignments = vec![None; n];
    let mut visited = vec![false; n];
    let mut cluster_id = 0;

    for i in 0..n {
        if visited[i] {
            continue;
        }
        visited[i] = true;

        // Find neighbors
        let neighbors: Vec<usize> = (0..n).filter(|&j| distances[[i, j]] <= eps).collect();

        if neighbors.len() < min_pts {
            // Mark as noise (will remain None)
            continue;
        }

        // Start new cluster
        let mut to_visit = neighbors.clone();
        cluster_assignments[i] = Some(cluster_id);

        let mut idx = 0;
        while idx < to_visit.len() {
            let point = to_visit[idx];

            if !visited[point] {
                visited[point] = true;

                let point_neighbors: Vec<usize> =
                    (0..n).filter(|&j| distances[[point, j]] <= eps).collect();

                if point_neighbors.len() >= min_pts {
                    // Add new neighbors to visit list
                    for &neighbor in &point_neighbors {
                        if !to_visit.contains(&neighbor) {
                            to_visit.push(neighbor);
                        }
                    }
                }
            }

            if cluster_assignments[point].is_none() {
                cluster_assignments[point] = Some(cluster_id);
            }

            idx += 1;
        }

        cluster_id += 1;
    }

    // Assign cluster IDs to minima
    for (i, minimum) in minima.iter_mut().enumerate() {
        minimum.cluster_id = cluster_assignments[i];
    }

    Ok(())
}

/// Simple threshold-based clustering
#[allow(dead_code)]
fn threshold_clustering<S>(
    minima: &mut [LocalMinimum<S>],
    options: &ClusteringOptions,
) -> Result<(), OptimizeError>
where
    S: Clone,
{
    let n = minima.len();
    if n <= 1 {
        if n == 1 {
            minima[0].cluster_id = Some(0);
        }
        return Ok(());
    }

    let mut cluster_assignments = vec![None; n];
    let mut cluster_id = 0;

    for i in 0..n {
        if cluster_assignments[i].is_some() {
            continue;
        }

        // Start new cluster
        cluster_assignments[i] = Some(cluster_id);

        // Find all points within threshold
        for j in (i + 1)..n {
            if cluster_assignments[j].is_some() {
                continue;
            }

            let distance = compute_distance(&minima[i], &minima[j], options);
            if distance <= options.distance_threshold {
                cluster_assignments[j] = Some(cluster_id);
            }
        }

        cluster_id += 1;
    }

    // Assign cluster IDs to minima
    for (i, minimum) in minima.iter_mut().enumerate() {
        minimum.cluster_id = cluster_assignments[i];
    }

    Ok(())
}

/// Compute distance matrix between all minima
#[allow(dead_code)]
fn compute_distance_matrix<S>(
    minima: &[LocalMinimum<S>],
    options: &ClusteringOptions,
) -> Array2<f64>
where
    S: Clone,
{
    let n = minima.len();
    let mut distances = Array2::zeros((n, n));

    for i in 0..n {
        for j in (i + 1)..n {
            let dist = compute_distance(&minima[i], &minima[j], options);
            distances[[i, j]] = dist;
            distances[[j, i]] = dist;
        }
    }

    distances
}

/// Compute distance between two local minima
#[allow(dead_code)]
fn compute_distance<S>(
    min1: &LocalMinimum<S>,
    min2: &LocalMinimum<S>,
    options: &ClusteringOptions,
) -> f64
where
    S: Clone,
{
    // Coordinate distance
    let coord_dist = (&min1.x - &min2.x).mapv(|x| x.powi(2)).sum().sqrt();

    if !options.use_function_values {
        return coord_dist;
    }

    // Function value distance
    let func_dist = (min1.f - min2.f).abs();

    // Combined distance
    coord_dist + options.function_weight * func_dist
}

/// Compute distance between two clusters (single linkage)
#[allow(dead_code)]
fn compute_cluster_distance(
    cluster1: &[usize],
    cluster2: &[usize],
    distances: &Array2<f64>,
) -> f64 {
    let mut min_dist = f64::INFINITY;

    for &i in cluster1 {
        for &j in cluster2 {
            let dist = distances[[i, j]];
            if dist < min_dist {
                min_dist = dist;
            }
        }
    }

    min_dist
}

/// Extract feature vectors for clustering
#[allow(dead_code)]
fn extract_features<S>(minima: &[LocalMinimum<S>], options: &ClusteringOptions) -> Array2<f64>
where
    S: Clone,
{
    let n = minima.len();
    if n == 0 {
        return Array2::zeros((0, 0));
    }

    let coord_dim = minima[0].x.len();
    let func_dim = if options.use_function_values { 1 } else { 0 };
    let total_dim = coord_dim + func_dim;

    let mut features = Array2::zeros((n, total_dim));

    // Extract coordinates
    for (i, minimum) in minima.iter().enumerate() {
        for j in 0..coord_dim {
            features[[i, j]] = minimum.x[j];
        }
    }

    // Add function values if requested
    if options.use_function_values {
        for (i, minimum) in minima.iter().enumerate() {
            features[[i, coord_dim]] = minimum.f * options.function_weight;
        }
    }

    // Normalize if requested
    if options.normalize_coordinates {
        normalize_features(&mut features, coord_dim);
    }

    features
}

/// Normalize feature matrix
#[allow(dead_code)]
fn normalize_features(features: &mut Array2<f64>, coord_dim: usize) {
    let n = features.nrows();
    if n == 0 {
        return;
    }

    // Normalize coordinates only
    for j in 0..coord_dim {
        let col = features.column(j);
        let min_val = col.iter().fold(f64::INFINITY, |a, &b| f64::min(a, b));
        let max_val = col.iter().fold(f64::NEG_INFINITY, |a, &b| f64::max(a, b));

        if (max_val - min_val).abs() > 1e-10 {
            for i in 0..n {
                features[[i, j]] = (features[[i, j]] - min_val) / (max_val - min_val);
            }
        }
    }
}

/// Initialize centroids using k-means++ algorithm
#[allow(dead_code)]
fn initialize_centroids_plus_plus(features: &Array2<f64>, k: usize) -> Array2<f64> {
    let (n, dim) = features.dim();
    let mut centroids = Array2::zeros((k, dim));

    if n == 0 || k == 0 {
        return centroids;
    }

    // Choose first centroid randomly
    let first_idx = 0; // In practice, should be random
    centroids.row_mut(0).assign(&features.row(first_idx));

    // Choose remaining centroids
    for c in 1..k {
        let mut distances = vec![f64::INFINITY; n];

        // Compute distance to nearest centroid for each point
        for (i, distance) in distances.iter_mut().enumerate().take(n) {
            let point = features.row(i);
            for j in 0..c {
                let centroid = centroids.row(j);
                let dist = euclidean_distance(&point, &centroid);
                *distance = distance.min(dist);
            }
        }

        // Choose next centroid (should be weighted random, using max for simplicity)
        let next_idx = distances
            .iter()
            .enumerate()
            .max_by(|a, b| a.1.partial_cmp(b.1).expect("Operation failed"))
            .map(|(i, _)| i)
            .unwrap_or(0);

        centroids.row_mut(c).assign(&features.row(next_idx));
    }

    centroids
}

/// Compute Euclidean distance between two points
#[allow(dead_code)]
fn euclidean_distance(a: &ArrayView1<f64>, b: &ArrayView1<f64>) -> f64 {
    (a - b).mapv(|x| x.powi(2)).sum().sqrt()
}

/// Compute cluster centroids and statistics
#[allow(dead_code)]
fn compute_cluster_centroids<S>(
    minima: &[LocalMinimum<S>],
) -> Result<Vec<ClusterCentroid>, OptimizeError>
where
    S: Clone,
{
    if minima.is_empty() {
        return Ok(Vec::new());
    }

    // Group minima by cluster
    let mut clusters: HashMap<usize, Vec<&LocalMinimum<S>>> = HashMap::new();

    for minimum in minima {
        if let Some(cluster_id) = minimum.cluster_id {
            clusters.entry(cluster_id).or_default().push(minimum);
        }
    }

    let mut centroids = Vec::new();

    for (_, cluster_minima) in clusters {
        if cluster_minima.is_empty() {
            continue;
        }

        let dim = cluster_minima[0].x.len();
        let mut centroid_x = Array1::zeros(dim);
        let mut f_sum = 0.0;
        let mut f_min = f64::INFINITY;

        // Compute centroid coordinates and function statistics
        for minimum in &cluster_minima {
            centroid_x = &centroid_x + &minimum.x;
            f_sum += minimum.f;
            f_min = f_min.min(minimum.f);
        }

        let size = cluster_minima.len();
        centroid_x /= size as f64;
        let f_avg = f_sum / size as f64;

        // Compute cluster radius
        let mut max_radius = 0.0;
        for minimum in &cluster_minima {
            let dist = (&minimum.x - &centroid_x).mapv(|x| x.powi(2)).sum().sqrt();
            max_radius = f64::max(max_radius, dist);
        }

        centroids.push(ClusterCentroid {
            x: centroid_x,
            f_avg,
            f_min,
            size,
            radius: max_radius,
        });
    }

    Ok(centroids)
}

/// Compute within-cluster sum of squares
#[allow(dead_code)]
fn compute_wcss<S>(minima: &[LocalMinimum<S>], centroids: &[ClusterCentroid]) -> f64
where
    S: Clone,
{
    let mut wcss = 0.0;

    for minimum in minima {
        if let Some(cluster_id) = minimum.cluster_id {
            if cluster_id < centroids.len() {
                let centroid = &centroids[cluster_id];
                let dist = (&minimum.x - &centroid.x).mapv(|x| x.powi(2)).sum();
                wcss += dist;
            }
        }
    }

    wcss
}

/// Compute silhouette score for clustering quality
#[allow(dead_code)]
fn compute_silhouette_score<S>(minima: &[LocalMinimum<S>]) -> Option<f64>
where
    S: Clone,
{
    if minima.len() < 2 {
        return None;
    }

    let mut silhouette_sum = 0.0;
    let mut valid_points = 0;

    for (i, minimum) in minima.iter().enumerate() {
        if let Some(cluster_id) = minimum.cluster_id {
            // Compute intra-cluster distance
            let mut intra_sum = 0.0;
            let mut intra_count = 0;

            // Compute inter-cluster distances
            let mut min_inter = f64::INFINITY;
            let mut cluster_inter_sums: HashMap<usize, f64> = HashMap::new();
            let mut cluster_inter_counts: HashMap<usize, usize> = HashMap::new();

            for (j, other) in minima.iter().enumerate() {
                if i == j {
                    continue;
                }

                if let Some(other_cluster_id) = other.cluster_id {
                    let dist = (&minimum.x - &other.x).mapv(|x| x.powi(2)).sum().sqrt();

                    if other_cluster_id == cluster_id {
                        intra_sum += dist;
                        intra_count += 1;
                    } else {
                        *cluster_inter_sums.entry(other_cluster_id).or_insert(0.0) += dist;
                        *cluster_inter_counts.entry(other_cluster_id).or_insert(0) += 1;
                    }
                }
            }

            // Find minimum inter-cluster distance
            for (other_cluster_id, sum) in cluster_inter_sums {
                let count = cluster_inter_counts[&other_cluster_id];
                if count > 0 {
                    let avg_inter = sum / count as f64;
                    min_inter = min_inter.min(avg_inter);
                }
            }

            if intra_count > 0 && min_inter < f64::INFINITY {
                let a = intra_sum / intra_count as f64;
                let b = min_inter;
                let silhouette = (b - a) / f64::max(a, b);
                silhouette_sum += silhouette;
                valid_points += 1;
            }
        }
    }

    if valid_points > 0 {
        Some(silhouette_sum / valid_points as f64)
    } else {
        None
    }
}

/// Generate diverse starting points for multi-start optimization
#[allow(dead_code)]
pub fn generate_diverse_start_points(
    bounds: &[(f64, f64)],
    num_points: usize,
    strategy: StartPointStrategy,
) -> Vec<Array1<f64>> {
    match strategy {
        StartPointStrategy::Random => generate_random_points(bounds, num_points),
        StartPointStrategy::LatinHypercube => generate_latin_hypercube_points(bounds, num_points),
        StartPointStrategy::Grid => generate_grid_points(bounds, num_points),
        StartPointStrategy::Sobol => generate_sobol_points(bounds, num_points),
    }
}

/// Strategy for generating starting points
#[derive(Debug, Clone, Copy)]
pub enum StartPointStrategy {
    Random,
    LatinHypercube,
    Grid,
    Sobol,
}

/// Generate random starting points
#[allow(dead_code)]
fn generate_random_points(bounds: &[(f64, f64)], num_points: usize) -> Vec<Array1<f64>> {
    let dim = bounds.len();
    let mut points = Vec::new();

    for _ in 0..num_points {
        let mut point = Array1::zeros(dim);
        for (i, &(low, high)) in bounds.iter().enumerate() {
            // Simple pseudo-random (in practice, use proper RNG)
            let t = (i * num_points + points.len()) as f64 / (num_points * dim) as f64;
            let random_val = (t * 17.0).fract(); // Simple pseudo-random
            point[i] = low + random_val * (high - low);
        }
        points.push(point);
    }

    points
}

/// Generate Latin Hypercube sampling points
#[allow(dead_code)]
fn generate_latin_hypercube_points(bounds: &[(f64, f64)], num_points: usize) -> Vec<Array1<f64>> {
    let dim = bounds.len();
    let mut points = Vec::new();

    // Simple Latin Hypercube implementation
    for i in 0..num_points {
        let mut point = Array1::zeros(dim);
        for (j, &(low, high)) in bounds.iter().enumerate() {
            let segment = (i as f64 + 0.5) / num_points as f64; // Center of segment
            point[j] = low + segment * (high - low);
        }
        points.push(point);
    }

    points
}

/// Generate grid points
#[allow(dead_code)]
fn generate_grid_points(bounds: &[(f64, f64)], num_points: usize) -> Vec<Array1<f64>> {
    let dim = bounds.len();
    if dim == 0 {
        return Vec::new();
    }

    let points_per_dim = ((num_points as f64).powf(1.0 / dim as f64)).ceil() as usize;
    let mut _points = Vec::new();

    // Generate grid coordinates recursively
    fn generate_grid_recursive(
        bounds: &[(f64, f64)],
        points_per_dim: usize,
        current_point: &mut Array1<f64>,
        dim_idx: usize,
        points: &mut Vec<Array1<f64>>,
    ) {
        if dim_idx >= bounds.len() {
            points.push(current_point.clone());
            return;
        }

        let (low, high) = bounds[dim_idx];
        for i in 0..points_per_dim {
            let t = if points_per_dim == 1 {
                0.5
            } else {
                i as f64 / (points_per_dim - 1) as f64
            };
            current_point[dim_idx] = low + t * (high - low);
            generate_grid_recursive(bounds, points_per_dim, current_point, dim_idx + 1, points);
        }
    }

    let mut current_point = Array1::zeros(dim);
    generate_grid_recursive(bounds, points_per_dim, &mut current_point, 0, &mut _points);

    // Truncate to requested number of _points
    _points.truncate(num_points);
    _points
}

/// Generate Sobol sequence points (simplified)
#[allow(dead_code)]
fn generate_sobol_points(bounds: &[(f64, f64)], num_points: usize) -> Vec<Array1<f64>> {
    // Simplified Sobol sequence (in practice, use proper implementation)
    let dim = bounds.len();
    let mut points = Vec::new();

    for i in 0..num_points {
        let mut point = Array1::zeros(dim);
        for (j, &(low, high)) in bounds.iter().enumerate() {
            // Van der Corput sequence for dimension j
            let mut n = i + 1;
            let base = 2 + j; // Different base for each dimension
            let mut result = 0.0;
            let mut f = 1.0 / base as f64;

            while n > 0 {
                result += (n % base) as f64 * f;
                n /= base;
                f /= base as f64;
            }

            point[j] = low + result * (high - low);
        }
        points.push(point);
    }

    points
}

#[cfg(test)]
mod tests {
    use super::*;
    use approx::assert_abs_diff_eq;

    #[test]
    fn test_simple_clustering() {
        // Create some test minima
        let mut minima = vec![
            LocalMinimum {
                x: Array1::from_vec(vec![0.0, 0.0]),
                f: 1.0,
                fun_value: 1.0,
                nit: 10,
                func_evals: 20,
                success: true,
                start_point: Array1::from_vec(vec![1.0, 1.0]),
                cluster_id: None,
                cluster_distance: None,
            },
            LocalMinimum {
                x: Array1::from_vec(vec![0.1, 0.1]),
                f: 1.1,
                fun_value: 1.1,
                nit: 12,
                func_evals: 24,
                success: true,
                start_point: Array1::from_vec(vec![1.1, 1.1]),
                cluster_id: None,
                cluster_distance: None,
            },
            LocalMinimum {
                x: Array1::from_vec(vec![5.0, 5.0]),
                f: 2.0,
                fun_value: 2.0,
                nit: 15,
                func_evals: 30,
                success: true,
                start_point: Array1::from_vec(vec![5.5, 5.5]),
                cluster_id: None,
                cluster_distance: None,
            },
        ];

        let options = ClusteringOptions {
            distance_threshold: 1.0,
            algorithm: ClusteringAlgorithm::Threshold,
            ..Default::default()
        };

        threshold_clustering(&mut minima, &options).expect("Operation failed");

        // First two should be in same cluster, third in different cluster
        assert_eq!(minima[0].cluster_id, minima[1].cluster_id);
        assert_ne!(minima[0].cluster_id, minima[2].cluster_id);
    }

    #[test]
    fn test_distance_computation() {
        let min1 = LocalMinimum {
            x: Array1::from_vec(vec![0.0, 0.0]),
            f: 1.0,
            fun_value: 1.0,
            nit: 10,
            func_evals: 20,
            success: true,
            start_point: Array1::from_vec(vec![1.0, 1.0]),
            cluster_id: None,
            cluster_distance: None,
        };

        let min2 = LocalMinimum {
            x: Array1::from_vec(vec![3.0, 4.0]),
            f: 2.0,
            fun_value: 2.0,
            nit: 12,
            func_evals: 24,
            success: true,
            start_point: Array1::from_vec(vec![3.5, 4.5]),
            cluster_id: None,
            cluster_distance: None,
        };

        let options = ClusteringOptions::default();
        let distance = compute_distance(&min1, &min2, &options);

        // Euclidean distance is 5.0, plus function weight * |1.0 - 2.0| = 5.0 + 0.1 * 1.0 = 5.1
        assert_abs_diff_eq!(distance, 5.1, epsilon = 1e-10);
    }

    #[test]
    fn test_start_point_generation() {
        let bounds = vec![(0.0, 10.0), (-5.0, 5.0)];
        let num_points = 5;

        let random_points =
            generate_diverse_start_points(&bounds, num_points, StartPointStrategy::Random);
        assert_eq!(random_points.len(), num_points);

        for point in &random_points {
            assert!(point[0] >= 0.0 && point[0] <= 10.0);
            assert!(point[1] >= -5.0 && point[1] <= 5.0);
        }

        let grid_points =
            generate_diverse_start_points(&bounds, num_points, StartPointStrategy::Grid);
        assert_eq!(grid_points.len(), num_points);

        for point in &grid_points {
            assert!(point[0] >= 0.0 && point[0] <= 10.0);
            assert!(point[1] >= -5.0 && point[1] <= 5.0);
        }
    }
}