scirs2_graph/generators/random_graphs.rs
1//! Advanced random graph models and network generation
2//!
3//! This module implements a suite of state-of-the-art random graph models:
4//!
5//! - **Erdős–Rényi G(n,p)**: edges exist independently with probability p
6//! - **Erdős–Rényi G(n,m)**: exactly m edges chosen uniformly at random
7//! - **Barabási–Albert**: preferential attachment produces scale-free networks
8//! - **Watts–Strogatz**: ring-lattice rewiring produces small-world topology
9//! - **Random d-regular**: uniform random d-regular graph via configuration model
10//! - **Hyperbolic random graph (HRG)**: geometric model in hyperbolic disk
11//! - **Stochastic Kronecker graph**: iterative tensor-product graph model
12//! - **Chung–Lu**: random graph with prescribed expected degree sequence
13
14use crate::base::Graph;
15use crate::error::{GraphError, Result};
16use scirs2_core::ndarray::Array2;
17use scirs2_core::rand_prelude::IndexedRandom;
18use scirs2_core::random::prelude::*;
19use scirs2_core::random::seq::SliceRandom;
20use std::collections::HashSet;
21
22// ─────────────────────────────────────────────────────────────────────────────
23// G(n, p) — Erdős–Rényi probability model
24// ─────────────────────────────────────────────────────────────────────────────
25
26/// Generate an Erdős–Rényi G(n, p) random graph.
27///
28/// Each pair of nodes (u, v) with u < v is connected independently with
29/// probability `p`. The expected number of edges is p · C(n, 2).
30///
31/// # Arguments
32/// * `n` – number of nodes (nodes are labelled 0 … n-1)
33/// * `p` – edge probability in \[0, 1\]
34/// * `rng` – a seeded or default random-number generator
35///
36/// # Errors
37/// Returns `GraphError::InvalidGraph` when `p ∉ [0,1]` or `n == 0`.
38///
39/// # Example
40/// ```rust
41/// use scirs2_graph::generators::random_graphs::erdos_renyi_g_np;
42/// use scirs2_core::random::prelude::*;
43/// let mut rng = StdRng::seed_from_u64(7);
44/// let g = erdos_renyi_g_np(20, 0.3, &mut rng).unwrap();
45/// assert_eq!(g.node_count(), 20);
46/// ```
47pub fn erdos_renyi_g_np<R: Rng>(n: usize, p: f64, rng: &mut R) -> Result<Graph<usize, f64>> {
48 if !(0.0..=1.0).contains(&p) {
49 return Err(GraphError::InvalidGraph(format!(
50 "erdos_renyi_g_np: p={p} must be in [0,1]"
51 )));
52 }
53 let mut g = Graph::new();
54 for i in 0..n {
55 g.add_node(i);
56 }
57 for u in 0..n {
58 for v in (u + 1)..n {
59 if rng.random::<f64>() < p {
60 g.add_edge(u, v, 1.0)?;
61 }
62 }
63 }
64 Ok(g)
65}
66
67// ─────────────────────────────────────────────────────────────────────────────
68// G(n, m) — Erdős–Rényi exact-edge-count model
69// ─────────────────────────────────────────────────────────────────────────────
70
71/// Generate an Erdős–Rényi G(n, m) random graph with **exactly** `m` edges.
72///
73/// A uniformly random subset of size `m` is chosen from all C(n, 2) possible
74/// edges using reservoir sampling (Fisher–Yates on the candidate list).
75///
76/// # Arguments
77/// * `n` – number of nodes
78/// * `m` – exact number of edges; must satisfy m ≤ C(n, 2)
79/// * `rng` – random-number generator
80///
81/// # Errors
82/// Returns `GraphError::InvalidGraph` when the edge count is infeasible.
83pub fn erdos_renyi_g_nm<R: Rng>(n: usize, m: usize, rng: &mut R) -> Result<Graph<usize, f64>> {
84 let max_edges = n.saturating_mul(n.saturating_sub(1)) / 2;
85 if m > max_edges {
86 return Err(GraphError::InvalidGraph(format!(
87 "erdos_renyi_g_nm: m={m} > C({n},2)={max_edges}"
88 )));
89 }
90 let mut g = Graph::new();
91 for i in 0..n {
92 g.add_node(i);
93 }
94 if m == 0 {
95 return Ok(g);
96 }
97
98 // Build the full candidate list and shuffle the first `m` entries.
99 let mut candidates: Vec<(usize, usize)> = Vec::with_capacity(max_edges);
100 for u in 0..n {
101 for v in (u + 1)..n {
102 candidates.push((u, v));
103 }
104 }
105 // Partial Fisher–Yates: select m edges
106 for i in 0..m {
107 let j = i + rng.random_range(0..(candidates.len() - i));
108 candidates.swap(i, j);
109 }
110 for &(u, v) in &candidates[..m] {
111 g.add_edge(u, v, 1.0)?;
112 }
113 Ok(g)
114}
115
116// ─────────────────────────────────────────────────────────────────────────────
117// Barabási–Albert preferential attachment
118// ─────────────────────────────────────────────────────────────────────────────
119
120/// Generate a Barabási–Albert (BA) scale-free graph via preferential attachment.
121///
122/// The model starts from an initial clique of `m + 1` nodes and then adds one
123/// new node at a time, connecting it to `m` existing nodes with probability
124/// proportional to their current degree (linear preferential attachment).
125///
126/// The resulting degree distribution follows a power law P(k) ~ k^{-3}.
127///
128/// # Arguments
129/// * `n` – total number of nodes (must satisfy n > m ≥ 1)
130/// * `m` – edges added per new node
131/// * `rng` – random-number generator
132///
133/// # Errors
134/// Returns `GraphError::InvalidGraph` for invalid parameter combinations.
135pub fn barabasi_albert<R: Rng>(n: usize, m: usize, rng: &mut R) -> Result<Graph<usize, f64>> {
136 if m == 0 {
137 return Err(GraphError::InvalidGraph(
138 "barabasi_albert: m must be ≥ 1".to_string(),
139 ));
140 }
141 if n <= m {
142 return Err(GraphError::InvalidGraph(format!(
143 "barabasi_albert: n={n} must be > m={m}"
144 )));
145 }
146
147 let mut g = Graph::new();
148
149 // Seed: complete graph on m+1 nodes
150 for i in 0..=m {
151 g.add_node(i);
152 }
153 for u in 0..=m {
154 for v in (u + 1)..=m {
155 g.add_edge(u, v, 1.0)?;
156 }
157 }
158
159 // Degree bookkeeping for O(1) preferential sampling via repeated-stub trick.
160 // stubs[i] appears degree[i] times so that uniform sampling over stubs gives
161 // the correct preferential-attachment distribution.
162 let mut stubs: Vec<usize> = Vec::with_capacity(n * m * 2);
163 for i in 0..=m {
164 for _ in 0..m {
165 stubs.push(i);
166 }
167 }
168
169 for new_node in (m + 1)..n {
170 g.add_node(new_node);
171 let mut chosen: HashSet<usize> = HashSet::with_capacity(m);
172
173 while chosen.len() < m {
174 let idx = rng.random_range(0..stubs.len());
175 let target = stubs[idx];
176 if target != new_node {
177 chosen.insert(target);
178 }
179 }
180
181 for &t in &chosen {
182 g.add_edge(new_node, t, 1.0)?;
183 // Update stubs: each accepted edge adds one stub for the target and one
184 // for the new node, preserving proportional-to-degree sampling.
185 stubs.push(t);
186 stubs.push(new_node);
187 }
188 }
189 Ok(g)
190}
191
192// ─────────────────────────────────────────────────────────────────────────────
193// Watts–Strogatz small-world model
194// ─────────────────────────────────────────────────────────────────────────────
195
196/// Generate a Watts–Strogatz small-world graph.
197///
198/// Starting from a ring lattice where every node is connected to its `k`
199/// nearest neighbours (k/2 on each side, k must be even), each edge is
200/// independently rewired with probability `beta` to a uniformly random
201/// non-neighbour. This interpolates between a regular lattice (β=0) and an
202/// Erdős–Rényi random graph (β=1) while passing through a small-world regime.
203///
204/// # Arguments
205/// * `n` – number of nodes
206/// * `k` – mean degree in the initial lattice (must be even, k < n)
207/// * `beta` – rewiring probability in \[0, 1\]
208/// * `rng` – random-number generator
209pub fn watts_strogatz<R: Rng>(
210 n: usize,
211 k: usize,
212 beta: f64,
213 rng: &mut R,
214) -> Result<Graph<usize, f64>> {
215 if k == 0 || k >= n || !k.is_multiple_of(2) {
216 return Err(GraphError::InvalidGraph(format!(
217 "watts_strogatz: k={k} must be even, ≥ 2, and < n={n}"
218 )));
219 }
220 if !(0.0..=1.0).contains(&beta) {
221 return Err(GraphError::InvalidGraph(format!(
222 "watts_strogatz: beta={beta} must be in [0,1]"
223 )));
224 }
225
226 // Adjacency set for O(1) duplicate-edge checking
227 let mut adj: Vec<HashSet<usize>> = vec![HashSet::new(); n];
228
229 // Initial ring lattice
230 for u in 0..n {
231 for s in 1..=(k / 2) {
232 let v = (u + s) % n;
233 adj[u].insert(v);
234 adj[v].insert(u);
235 }
236 }
237
238 // Rewiring pass: iterate over u, then over each right-half neighbour v.
239 for u in 0..n {
240 for s in 1..=(k / 2) {
241 let v = (u + s) % n;
242 if !adj[u].contains(&v) {
243 // Edge was already rewired away in a previous step; skip
244 continue;
245 }
246 if rng.random::<f64>() < beta {
247 // Choose a new target w != u and not already a neighbour
248 let mut w = rng.random_range(0..n);
249 let mut attempts = 0usize;
250 while (w == u || adj[u].contains(&w)) && attempts < n * 4 {
251 w = rng.random_range(0..n);
252 attempts += 1;
253 }
254 if w == u || adj[u].contains(&w) {
255 // Could not find a valid rewire target; skip
256 continue;
257 }
258 // Remove edge (u, v)
259 adj[u].remove(&v);
260 adj[v].remove(&u);
261 // Add edge (u, w)
262 adj[u].insert(w);
263 adj[w].insert(u);
264 }
265 }
266 }
267
268 // Build the Graph from the final adjacency sets
269 let mut g = Graph::new();
270 for i in 0..n {
271 g.add_node(i);
272 }
273 for u in 0..n {
274 for &v in &adj[u] {
275 if u < v {
276 g.add_edge(u, v, 1.0)?;
277 }
278 }
279 }
280
281 Ok(g)
282}
283
284// ─────────────────────────────────────────────────────────────────────────────
285// Random d-regular graph
286// ─────────────────────────────────────────────────────────────────────────────
287
288/// Generate a uniformly random d-regular graph on `n` nodes.
289///
290/// A *d-regular* graph is one where every node has exactly degree `d`.
291/// This is feasible only when `n·d` is even.
292///
293/// The algorithm uses the **configuration model** with self-loop / parallel-edge
294/// rejection: stubs are generated (n·d stubs total), randomly paired, and
295/// the pairing is rejected and retried if it produces a self-loop or
296/// multi-edge. The expected number of retries is O(1) for fixed d.
297///
298/// # Arguments
299/// * `n` – number of nodes
300/// * `d` – required degree of every node
301/// * `rng` – random-number generator
302///
303/// # Returns
304/// `Some(Graph)` when a valid d-regular graph was found within the attempt
305/// budget, `None` when the degree sequence is infeasible or sampling
306/// consistently fails (e.g., n < d+1).
307pub fn random_regular<R: Rng>(n: usize, d: usize, rng: &mut R) -> Option<Graph<usize, f64>> {
308 if n == 0 || d == 0 {
309 let mut g = Graph::new();
310 for i in 0..n {
311 g.add_node(i);
312 }
313 return Some(g);
314 }
315 if d >= n {
316 return None; // Not possible without self-loops
317 }
318 if !(n * d).is_multiple_of(2) {
319 return None; // Degree sequence not graphical
320 }
321
322 let max_outer_attempts = 100usize;
323
324 for _ in 0..max_outer_attempts {
325 if let Some(g) = try_build_regular(n, d, rng) {
326 return Some(g);
327 }
328 }
329 None
330}
331
332/// Single attempt to build a d-regular graph via the configuration model.
333fn try_build_regular<R: Rng>(n: usize, d: usize, rng: &mut R) -> Option<Graph<usize, f64>> {
334 // Create stubs: node i appears d times
335 let mut stubs: Vec<usize> = (0..n).flat_map(|i| std::iter::repeat_n(i, d)).collect();
336 stubs.shuffle(rng);
337
338 let mut adj: Vec<HashSet<usize>> = vec![HashSet::new(); n];
339
340 while stubs.len() >= 2 {
341 let u = stubs[0];
342 let v = stubs[1];
343 if u == v || adj[u].contains(&v) {
344 // Try to find a valid partner for u via random swap
345 let swap_idx = 2 + rng.random_range(0..stubs.len().saturating_sub(2).max(1));
346 if swap_idx < stubs.len() {
347 stubs.swap(1, swap_idx);
348 } else {
349 return None; // Give up on this attempt
350 }
351 // Attempt limit guard: after many swaps we give up
352 continue;
353 }
354 adj[u].insert(v);
355 adj[v].insert(u);
356 stubs.drain(0..2);
357 }
358
359 if !stubs.is_empty() {
360 return None;
361 }
362
363 let mut g = Graph::new();
364 for i in 0..n {
365 g.add_node(i);
366 }
367 for u in 0..n {
368 for &v in &adj[u] {
369 if u < v {
370 g.add_edge(u, v, 1.0).ok()?;
371 }
372 }
373 }
374 Some(g)
375}
376
377// ─────────────────────────────────────────────────────────────────────────────
378// Hyperbolic Random Graph (HRG)
379// ─────────────────────────────────────────────────────────────────────────────
380
381/// Generate a hyperbolic random graph (HRG) in the Poincaré disk model.
382///
383/// `n` nodes are placed uniformly at random in a hyperbolic disk of radius `r`
384/// according to the quasi-uniform distribution with curvature parameter `alpha`
385/// (controls degree heterogeneity; `alpha = 0.5` recovers the pure geometric
386/// model, higher values concentrate nodes near the boundary). Two nodes are
387/// connected if their hyperbolic distance is at most `r`.
388///
389/// The hyperbolic distance between nodes at (r₁, θ₁) and (r₂, θ₂) in
390/// polar-hyperbolic coordinates is:
391///
392/// ```text
393/// d(u,v) = acosh(cosh(r₁)·cosh(r₂) − sinh(r₁)·sinh(r₂)·cos(θ₁−θ₂))
394/// ```
395///
396/// # Arguments
397/// * `n` – number of nodes
398/// * `r` – disk radius; larger `r` gives sparser graphs
399/// * `alpha` – radial density exponent (> 0; 0.5 ≤ α ≤ 1 for scale-free degree)
400/// * `rng` – random-number generator
401///
402/// # Reference
403/// Krioukov et al., "Hyperbolic geometry of complex networks", Phys. Rev. E,
404/// 82(3), 036106, 2010.
405pub fn hyperbolic_random_graph<R: Rng>(
406 n: usize,
407 r: f64,
408 alpha: f64,
409 rng: &mut R,
410) -> Result<Graph<usize, f64>> {
411 if r <= 0.0 {
412 return Err(GraphError::InvalidGraph(format!(
413 "hyperbolic_random_graph: r={r} must be > 0"
414 )));
415 }
416 if alpha <= 0.0 {
417 return Err(GraphError::InvalidGraph(format!(
418 "hyperbolic_random_graph: alpha={alpha} must be > 0"
419 )));
420 }
421
422 let mut g = Graph::new();
423 for i in 0..n {
424 g.add_node(i);
425 }
426 if n < 2 {
427 return Ok(g);
428 }
429
430 // Sample polar coordinates in the hyperbolic disk.
431 // Radial coordinate: the quasi-uniform distribution on [0, R] is
432 // F(ρ) = (cosh(alpha·ρ) − 1) / (cosh(alpha·R) − 1)
433 // We invert it via inverse-CDF: ρ = acosh(1 + u·(cosh(α·R)−1)) / α
434 let cosh_alpha_r = (alpha * r).cosh();
435 let two_pi = std::f64::consts::TAU;
436
437 let mut coords: Vec<(f64, f64)> = Vec::with_capacity(n); // (rho, theta)
438 for _ in 0..n {
439 let u: f64 = rng.random();
440 let rho = ((1.0 + u * (cosh_alpha_r - 1.0)).acosh()) / alpha;
441 let theta = rng.random::<f64>() * two_pi;
442 coords.push((rho, theta));
443 }
444
445 // Connect pairs within hyperbolic distance r
446 for i in 0..n {
447 let (r1, t1) = coords[i];
448 let (sinh_r1, cosh_r1) = (r1.sinh(), r1.cosh());
449 for j in (i + 1)..n {
450 let (r2, t2) = coords[j];
451 let delta_theta = (t1 - t2).abs();
452 // Map delta_theta to [0, π]
453 let delta_theta = if delta_theta > std::f64::consts::PI {
454 two_pi - delta_theta
455 } else {
456 delta_theta
457 };
458 let arg = cosh_r1 * r2.cosh() - sinh_r1 * r2.sinh() * delta_theta.cos();
459 // arg can be < 1 due to floating-point; clamp to avoid NaN from acosh
460 let hyp_dist = if arg <= 1.0 { 0.0 } else { arg.acosh() };
461 if hyp_dist <= r {
462 g.add_edge(i, j, hyp_dist)?;
463 }
464 }
465 }
466 Ok(g)
467}
468
469// ─────────────────────────────────────────────────────────────────────────────
470// Stochastic Kronecker graph
471// ─────────────────────────────────────────────────────────────────────────────
472
473/// Generate a stochastic Kronecker graph.
474///
475/// The *Kronecker graph* model (Leskovec et al., 2010) starts from a small
476/// `k₀ × k₀` *initiator matrix* Θ whose entries are probabilities in (0,1].
477/// The Kronecker product Θ^{⊗k} gives a (k₀^k × k₀^k) matrix of edge
478/// probabilities; an edge (i, j) is present independently with that probability.
479///
480/// This function computes the Kronecker product iteratively and samples edges
481/// without materialising the full N×N matrix: for each entry (i, j) of the
482/// product matrix, the entry value equals the product of initiator entries
483/// indexed by the base-k₀ digits of i and j.
484///
485/// # Arguments
486/// * `initiator` – square probability matrix of shape k₀ × k₀ (entries in (0, 1])
487/// * `k` – number of Kronecker iterations; produces N = k₀^k nodes
488/// * `rng` – random-number generator
489///
490/// # Errors
491/// Returns `GraphError::InvalidGraph` if the initiator is not square or contains
492/// values outside (0, 1].
493///
494/// # Reference
495/// Leskovec, J., Chakrabarti, D., Kleinberg, J., Faloutsos, C., & Ghahramani, Z.
496/// "Kronecker graphs: An approach to modeling networks." JMLR 11 (2010).
497pub fn kronecker_graph<R: Rng>(
498 initiator: &Array2<f64>,
499 k: usize,
500 rng: &mut R,
501) -> Result<Graph<usize, f64>> {
502 let k0 = initiator.nrows();
503 if k0 == 0 || initiator.ncols() != k0 {
504 return Err(GraphError::InvalidGraph(
505 "kronecker_graph: initiator must be a non-empty square matrix".to_string(),
506 ));
507 }
508 for val in initiator.iter() {
509 if !(*val > 0.0 && *val <= 1.0) {
510 return Err(GraphError::InvalidGraph(format!(
511 "kronecker_graph: initiator entries must be in (0,1], found {val}"
512 )));
513 }
514 }
515 if k == 0 {
516 return Err(GraphError::InvalidGraph(
517 "kronecker_graph: k must be ≥ 1".to_string(),
518 ));
519 }
520
521 let n = k0.pow(k as u32);
522 let mut g = Graph::new();
523 for i in 0..n {
524 g.add_node(i);
525 }
526
527 // Compute edge probability for (i,j) by decomposing i,j in base k0
528 for i in 0..n {
529 for j in (i + 1)..n {
530 let prob = kronecker_edge_prob(initiator, k0, k, i, j);
531 if rng.random::<f64>() < prob {
532 g.add_edge(i, j, 1.0)?;
533 }
534 }
535 }
536 Ok(g)
537}
538
539/// Compute the edge probability P(i,j) in the k-th Kronecker power.
540///
541/// The probability equals the product of initiator\[i_t\]\[j_t\] where
542/// (i_0, …, i_{k-1}) and (j_0, …, j_{k-1}) are the base-k₀ representations
543/// of i and j (most-significant digit first).
544fn kronecker_edge_prob(
545 initiator: &Array2<f64>,
546 k0: usize,
547 k: usize,
548 mut i: usize,
549 mut j: usize,
550) -> f64 {
551 let mut prob = 1.0f64;
552 for _ in 0..k {
553 let digit_i = i % k0;
554 let digit_j = j % k0;
555 prob *= initiator[[digit_i, digit_j]];
556 i /= k0;
557 j /= k0;
558 }
559 prob
560}
561
562// ─────────────────────────────────────────────────────────────────────────────
563// Chung–Lu model
564// ─────────────────────────────────────────────────────────────────────────────
565
566/// Generate a Chung–Lu random graph with prescribed expected degree sequence.
567///
568/// In the Chung–Lu model each node `i` has a weight `w_i` (intended expected
569/// degree). The probability of an edge (i, j) is:
570///
571/// ```text
572/// P(i,j) = min(w_i · w_j / W, 1) where W = Σ w_i
573/// ```
574///
575/// When the weights equal the observed degrees of a real network, this model
576/// reproduces the degree sequence in expectation and is often used as a
577/// null-model baseline.
578///
579/// # Arguments
580/// * `degree_seq` – slice of non-negative expected degrees; all values ≥ 0
581/// * `rng` – random-number generator
582///
583/// # Errors
584/// Returns `GraphError::InvalidGraph` if the degree sequence is empty or any
585/// value is negative.
586///
587/// # Reference
588/// Chung, F., & Lu, L. "Connected components in random graphs with given expected
589/// degree sequences." Annals of Combinatorics, 6(2), 125–145, 2002.
590pub fn chung_lu<R: Rng>(degree_seq: &[f64], rng: &mut R) -> Result<Graph<usize, f64>> {
591 if degree_seq.is_empty() {
592 return Ok(Graph::new());
593 }
594 for &w in degree_seq {
595 if w < 0.0 || !w.is_finite() {
596 return Err(GraphError::InvalidGraph(format!(
597 "chung_lu: degree sequence entries must be non-negative finite, found {w}"
598 )));
599 }
600 }
601
602 let n = degree_seq.len();
603 let total_weight: f64 = degree_seq.iter().sum();
604
605 let mut g = Graph::new();
606 for i in 0..n {
607 g.add_node(i);
608 }
609
610 if total_weight <= 0.0 {
611 return Ok(g);
612 }
613
614 for i in 0..n {
615 for j in (i + 1)..n {
616 let p = (degree_seq[i] * degree_seq[j] / total_weight).min(1.0);
617 if p > 0.0 && rng.random::<f64>() < p {
618 g.add_edge(i, j, 1.0)?;
619 }
620 }
621 }
622 Ok(g)
623}
624
625// ─────────────────────────────────────────────────────────────────────────────
626// Tests
627// ─────────────────────────────────────────────────────────────────────────────
628
629#[cfg(test)]
630mod tests {
631 use super::*;
632 use scirs2_core::random::prelude::StdRng;
633 use scirs2_core::random::SeedableRng;
634
635 // ── G(n,p) ──────────────────────────────────────────────────────────────
636
637 #[test]
638 fn test_gnp_basic() {
639 let mut rng = StdRng::seed_from_u64(1);
640 let g = erdos_renyi_g_np(10, 0.5, &mut rng).expect("gnp failed");
641 assert_eq!(g.node_count(), 10);
642 assert!(g.edge_count() <= 45);
643 }
644
645 #[test]
646 fn test_gnp_p0_no_edges() {
647 let mut rng = StdRng::seed_from_u64(2);
648 let g = erdos_renyi_g_np(15, 0.0, &mut rng).expect("gnp p=0 failed");
649 assert_eq!(g.edge_count(), 0);
650 }
651
652 #[test]
653 fn test_gnp_p1_complete() {
654 let mut rng = StdRng::seed_from_u64(3);
655 let g = erdos_renyi_g_np(8, 1.0, &mut rng).expect("gnp p=1 failed");
656 assert_eq!(g.edge_count(), 8 * 7 / 2);
657 }
658
659 #[test]
660 fn test_gnp_invalid_p() {
661 let mut rng = StdRng::seed_from_u64(4);
662 assert!(erdos_renyi_g_np(10, 1.5, &mut rng).is_err());
663 assert!(erdos_renyi_g_np(10, -0.1, &mut rng).is_err());
664 }
665
666 // ── G(n,m) ──────────────────────────────────────────────────────────────
667
668 #[test]
669 fn test_gnm_exact_edges() {
670 let mut rng = StdRng::seed_from_u64(5);
671 for m in [0, 1, 10, 45] {
672 let g = erdos_renyi_g_nm(10, m, &mut rng).expect("gnm failed");
673 assert_eq!(g.node_count(), 10);
674 assert_eq!(g.edge_count(), m);
675 }
676 }
677
678 #[test]
679 fn test_gnm_too_many_edges() {
680 let mut rng = StdRng::seed_from_u64(6);
681 assert!(erdos_renyi_g_nm(5, 11, &mut rng).is_err()); // C(5,2) = 10
682 }
683
684 // ── Barabási–Albert ──────────────────────────────────────────────────────
685
686 #[test]
687 fn test_ba_node_count() {
688 let mut rng = StdRng::seed_from_u64(7);
689 let g = barabasi_albert(50, 3, &mut rng).expect("ba failed");
690 assert_eq!(g.node_count(), 50);
691 }
692
693 #[test]
694 fn test_ba_min_edges() {
695 let mut rng = StdRng::seed_from_u64(8);
696 // Initial clique: C(m+1, 2) edges; each new node adds at least m edges.
697 let g = barabasi_albert(20, 2, &mut rng).expect("ba failed");
698 // seed clique has C(3,2)=3 edges + 17 * 2 = 37 total minimum
699 assert!(g.edge_count() >= 3 + 17 * 2);
700 }
701
702 #[test]
703 fn test_ba_invalid_params() {
704 let mut rng = StdRng::seed_from_u64(9);
705 assert!(barabasi_albert(5, 0, &mut rng).is_err());
706 assert!(barabasi_albert(3, 5, &mut rng).is_err());
707 }
708
709 // ── Watts–Strogatz ───────────────────────────────────────────────────────
710
711 #[test]
712 fn test_ws_lattice() {
713 let mut rng = StdRng::seed_from_u64(10);
714 let g = watts_strogatz(20, 4, 0.0, &mut rng).expect("ws failed");
715 assert_eq!(g.node_count(), 20);
716 // Each node connects to k=4 others in a ring, so n*k/2 edges
717 assert_eq!(g.edge_count(), 20 * 4 / 2);
718 }
719
720 #[test]
721 fn test_ws_random() {
722 let mut rng = StdRng::seed_from_u64(11);
723 let g = watts_strogatz(30, 4, 1.0, &mut rng).expect("ws full rewire failed");
724 assert_eq!(g.node_count(), 30);
725 // rewiring preserves number of edges (at most n*k/2)
726 assert!(g.edge_count() <= 30 * 4 / 2);
727 }
728
729 #[test]
730 fn test_ws_invalid_params() {
731 let mut rng = StdRng::seed_from_u64(12);
732 assert!(watts_strogatz(10, 3, 0.5, &mut rng).is_err()); // k odd
733 assert!(watts_strogatz(10, 10, 0.5, &mut rng).is_err()); // k >= n
734 assert!(watts_strogatz(10, 4, -0.1, &mut rng).is_err());
735 }
736
737 // ── Random regular ───────────────────────────────────────────────────────
738
739 #[test]
740 fn test_random_regular_degrees() {
741 let mut rng = StdRng::seed_from_u64(13);
742 if let Some(g) = random_regular(10, 3, &mut rng) {
743 assert_eq!(g.node_count(), 10);
744 // Every node should have degree 3
745 for i in 0..10usize {
746 assert_eq!(g.degree(&i), 3);
747 }
748 }
749 // n·d must be even: 10*3=30 is even, so it should succeed
750 }
751
752 #[test]
753 fn test_random_regular_infeasible() {
754 let mut rng = StdRng::seed_from_u64(14);
755 // n*d odd → infeasible
756 assert!(random_regular(5, 3, &mut rng).is_none()); // 15 is odd
757 // d >= n → infeasible
758 assert!(random_regular(4, 4, &mut rng).is_none());
759 }
760
761 // ── Hyperbolic random graph ───────────────────────────────────────────────
762
763 #[test]
764 fn test_hrg_node_count() {
765 let mut rng = StdRng::seed_from_u64(15);
766 let g = hyperbolic_random_graph(30, 6.0, 0.75, &mut rng).expect("hrg failed");
767 assert_eq!(g.node_count(), 30);
768 }
769
770 #[test]
771 fn test_hrg_small_radius_sparse() {
772 let mut rng = StdRng::seed_from_u64(16);
773 // Very small disk radius → very few edges (most pairs far apart)
774 let g = hyperbolic_random_graph(50, 0.5, 0.75, &mut rng).expect("hrg small r failed");
775 assert!(g.edge_count() < 50 * 49 / 2);
776 }
777
778 #[test]
779 fn test_hrg_invalid_params() {
780 let mut rng = StdRng::seed_from_u64(17);
781 assert!(hyperbolic_random_graph(10, -1.0, 0.5, &mut rng).is_err());
782 assert!(hyperbolic_random_graph(10, 5.0, 0.0, &mut rng).is_err());
783 }
784
785 // ── Stochastic Kronecker ─────────────────────────────────────────────────
786
787 #[test]
788 fn test_kronecker_small() {
789 use scirs2_core::ndarray::array;
790 let mut rng = StdRng::seed_from_u64(18);
791 // 2×2 initiator, k=3 → 2^3=8 nodes
792 let theta = array![[0.9, 0.5], [0.5, 0.3]];
793 let g = kronecker_graph(&theta, 3, &mut rng).expect("kronecker failed");
794 assert_eq!(g.node_count(), 8);
795 // All edge weights should be 1.0 (binary)
796 for e in g.edges() {
797 assert!((e.weight - 1.0).abs() < 1e-9);
798 }
799 }
800
801 #[test]
802 fn test_kronecker_invalid() {
803 use scirs2_core::ndarray::array;
804 let mut rng = StdRng::seed_from_u64(19);
805 // Non-square initiator
806 let non_square = Array2::zeros((2, 3));
807 assert!(kronecker_graph(&non_square, 2, &mut rng).is_err());
808 // k=0
809 let theta = array![[0.5, 0.5], [0.5, 0.5]];
810 assert!(kronecker_graph(&theta, 0, &mut rng).is_err());
811 }
812
813 // ── Chung–Lu ─────────────────────────────────────────────────────────────
814
815 #[test]
816 fn test_chung_lu_basic() {
817 let mut rng = StdRng::seed_from_u64(20);
818 let w = vec![5.0, 4.0, 3.0, 3.0, 2.0, 2.0, 1.0, 1.0];
819 let g = chung_lu(&w, &mut rng).expect("chung_lu failed");
820 assert_eq!(g.node_count(), 8);
821 // edges should exist (high weights ≈ high connectivity)
822 assert!(g.edge_count() > 0);
823 }
824
825 #[test]
826 fn test_chung_lu_zero_weights() {
827 let mut rng = StdRng::seed_from_u64(21);
828 let w = vec![0.0; 10];
829 let g = chung_lu(&w, &mut rng).expect("chung_lu zero weights failed");
830 assert_eq!(g.node_count(), 10);
831 assert_eq!(g.edge_count(), 0);
832 }
833
834 #[test]
835 fn test_chung_lu_invalid() {
836 let mut rng = StdRng::seed_from_u64(22);
837 assert!(chung_lu(&[-1.0, 2.0], &mut rng).is_err());
838 assert!(chung_lu(&[f64::INFINITY, 1.0], &mut rng).is_err());
839 }
840
841 #[test]
842 fn test_chung_lu_empty() {
843 let mut rng = StdRng::seed_from_u64(23);
844 let g = chung_lu(&[], &mut rng).expect("chung_lu empty failed");
845 assert_eq!(g.node_count(), 0);
846 }
847}