scirs2-io 0.5.0

Input/Output utilities module for SciRS2 (scirs2-io)
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
//! Matrix Market file format example
//!
//! This example demonstrates the Matrix Market file format capabilities:
//! - Creating sparse matrices in Matrix Market format
//! - Writing and reading Matrix Market files
//! - Working with different matrix types and symmetries
//! - Converting between different representations

use scirs2_core::ndarray::{Array1, Array2};
use scirs2_io::error::Result;
use scirs2_io::matrix_market::{
    coo_to_sparse, read_dense_matrix, read_sparse_matrix, sparse_to_coo, write_dense_matrix,
    write_sparse_matrix, MMDataType, MMDenseMatrix, MMFormat, MMHeader, MMSparseMatrix, MMSymmetry,
    SparseEntry,
};
use std::fs;

#[allow(dead_code)]
fn main() -> Result<()> {
    println!("=== Matrix Market Example ===");

    // Example 1: Create and write a sparse matrix
    create_sparse_matrix_example()?;

    // Example 2: Create and write a dense matrix
    create_dense_matrix_example()?;

    // Example 3: Read and analyze matrices
    read_and_analyze_matrices()?;

    // Example 4: Demonstrate different matrix types
    demonstrate_matrix_types()?;

    // Example 5: Coordinate format conversions
    demonstrate_coordinate_conversions()?;

    println!("Matrix Market example completed successfully!");
    Ok(())
}

#[allow(dead_code)]
fn create_sparse_matrix_example() -> Result<()> {
    println!("\n1. Creating sparse matrix example...");

    // Create a sparse matrix header
    let header = MMHeader {
        object: "matrix".to_string(),
        format: MMFormat::Coordinate,
        data_type: MMDataType::Real,
        symmetry: MMSymmetry::General,
        comments: vec![
            "Example sparse matrix".to_string(),
            "Generated by scirs2-io Matrix Market example".to_string(),
        ],
    };

    // Create a 5x5 sparse matrix with some non-zero entries
    let mut entries = Vec::new();

    // Diagonal entries
    for i in 0..5 {
        entries.push(SparseEntry {
            row: i,
            col: i,
            value: (i + 1) as f64,
        });
    }

    // Some off-diagonal entries
    entries.push(SparseEntry {
        row: 0,
        col: 2,
        value: 0.5,
    });
    entries.push(SparseEntry {
        row: 1,
        col: 3,
        value: -1.2,
    });
    entries.push(SparseEntry {
        row: 2,
        col: 4,
        value: 2.3,
    });
    entries.push(SparseEntry {
        row: 3,
        col: 0,
        value: -0.8,
    });
    entries.push(SparseEntry {
        row: 4,
        col: 1,
        value: 1.7,
    });

    let sparse_matrix = MMSparseMatrix {
        header,
        rows: 5,
        cols: 5,
        nnz: entries.len(),
        entries,
    };

    // Write the sparse matrix to a file
    write_sparse_matrix("example_sparse.mtx", &sparse_matrix)?;

    println!("  Created sparse matrix:");
    println!("    Size: {}x{}", sparse_matrix.rows, sparse_matrix.cols);
    println!("    Non-zeros: {}", sparse_matrix.nnz);
    println!("    Format: {:?}", sparse_matrix.header.format);
    println!("    Data type: {:?}", sparse_matrix.header.data_type);
    println!("    Symmetry: {:?}", sparse_matrix.header.symmetry);
    println!("    Written to: example_sparse.mtx");

    Ok(())
}

#[allow(dead_code)]
fn create_dense_matrix_example() -> Result<()> {
    println!("\n2. Creating dense matrix example...");

    // Create a dense matrix header
    let header = MMHeader {
        object: "matrix".to_string(),
        format: MMFormat::Array,
        data_type: MMDataType::Real,
        symmetry: MMSymmetry::General,
        comments: vec![
            "Example dense matrix".to_string(),
            "3x3 test matrix with various values".to_string(),
        ],
    };

    // Create a 3x3 dense matrix
    #[rustfmt::skip]
    let data = Array2::from_shape_vec(
        (3, 3),
        vec![
            1.0, 2.0, 3.0,  // First row
            4.0, 5.0, 6.0,  // Second row
            7.0, 8.0, 9.0,  // Third row
        ]
    ).expect("Operation failed");

    let dense_matrix = MMDenseMatrix {
        header,
        rows: 3,
        cols: 3,
        data,
    };

    // Write the dense matrix to a file
    write_dense_matrix("example_dense.mtx", &dense_matrix)?;

    println!("  Created dense matrix:");
    println!("    Size: {}x{}", dense_matrix.rows, dense_matrix.cols);
    println!("    Format: {:?}", dense_matrix.header.format);
    println!("    Data type: {:?}", dense_matrix.header.data_type);
    println!("    Matrix:");
    for row in 0..dense_matrix.rows {
        print!("      [");
        for col in 0..dense_matrix.cols {
            print!("{:6.2}", dense_matrix.data[[row, col]]);
            if col < dense_matrix.cols - 1 {
                print!(", ");
            }
        }
        println!("]");
    }
    println!("    Written to: example_dense.mtx");

    Ok(())
}

#[allow(dead_code)]
fn read_and_analyze_matrices() -> Result<()> {
    println!("\n3. Reading and analyzing matrices...");

    // Read the sparse matrix
    println!("  Reading sparse matrix...");
    let sparse_matrix = read_sparse_matrix("example_sparse.mtx")?;

    println!("    Loaded sparse matrix:");
    println!("      Size: {}x{}", sparse_matrix.rows, sparse_matrix.cols);
    println!("      Non-zeros: {}", sparse_matrix.nnz);
    println!("      Comments: {:?}", sparse_matrix.header.comments);

    // Analyze sparsity
    let total_elements = sparse_matrix.rows * sparse_matrix.cols;
    let sparsity = 1.0 - (sparse_matrix.nnz as f64 / total_elements as f64);
    println!(
        "      Sparsity: {:.2}% ({} zeros out of {} total elements)",
        sparsity * 100.0,
        total_elements - sparse_matrix.nnz,
        total_elements
    );

    // Show some entries
    println!("      First 5 non-zero entries:");
    for entry in sparse_matrix.entries.iter().take(5) {
        println!(
            "        ({}, {}) = {}",
            entry.row + 1,
            entry.col + 1,
            entry.value
        );
    }

    // Read the dense matrix
    println!("  Reading dense matrix...");
    let dense_matrix = read_dense_matrix("example_dense.mtx")?;

    println!("    Loaded dense matrix:");
    println!("      Size: {}x{}", dense_matrix.rows, dense_matrix.cols);
    println!("      Comments: {:?}", dense_matrix.header.comments);
    println!("      Data shape: {:?}", dense_matrix.data.shape());

    // Calculate some statistics
    let sum = dense_matrix.data.sum();
    let mean = sum / (dense_matrix.rows * dense_matrix.cols) as f64;
    let min = dense_matrix
        .data
        .iter()
        .fold(f64::INFINITY, |a, &b| a.min(b));
    let max = dense_matrix
        .data
        .iter()
        .fold(f64::NEG_INFINITY, |a, &b| a.max(b));

    println!("      Statistics:");
    println!("        Sum: {:.2}", sum);
    println!("        Mean: {:.2}", mean);
    println!("        Min: {:.2}", min);
    println!("        Max: {:.2}", max);

    Ok(())
}

#[allow(dead_code)]
fn demonstrate_matrix_types() -> Result<()> {
    println!("\n4. Demonstrating different matrix types...");

    // Create a symmetric matrix
    println!("  Creating symmetric matrix...");
    let sym_header = MMHeader {
        object: "matrix".to_string(),
        format: MMFormat::Coordinate,
        data_type: MMDataType::Real,
        symmetry: MMSymmetry::Symmetric,
        comments: vec!["Symmetric matrix example".to_string()],
    };

    // For symmetric matrices, we only need to store the lower triangle
    let sym_entries = vec![
        SparseEntry {
            row: 0,
            col: 0,
            value: 2.0,
        }, // Diagonal
        SparseEntry {
            row: 1,
            col: 0,
            value: 1.0,
        }, // Lower triangle only
        SparseEntry {
            row: 1,
            col: 1,
            value: 3.0,
        }, // Diagonal
        SparseEntry {
            row: 2,
            col: 0,
            value: 0.5,
        }, // Lower triangle only
        SparseEntry {
            row: 2,
            col: 1,
            value: -0.8,
        }, // Lower triangle only
        SparseEntry {
            row: 2,
            col: 2,
            value: 4.0,
        }, // Diagonal
    ];

    let sym_matrix = MMSparseMatrix {
        header: sym_header,
        rows: 3,
        cols: 3,
        nnz: sym_entries.len(),
        entries: sym_entries,
    };

    write_sparse_matrix("symmetric_matrix.mtx", &sym_matrix)?;
    println!(
        "    Written symmetric matrix (3x3) with {} stored entries",
        sym_matrix.nnz
    );

    // Create a pattern matrix (no values, just structure)
    println!("  Creating pattern matrix...");
    let pattern_header = MMHeader {
        object: "matrix".to_string(),
        format: MMFormat::Coordinate,
        data_type: MMDataType::Pattern,
        symmetry: MMSymmetry::General,
        comments: vec!["Pattern matrix - structure only".to_string()],
    };

    // Pattern matrices don't store values, just positions
    let pattern_entries = vec![
        SparseEntry {
            row: 0,
            col: 0,
            value: 1.0,
        }, // Value ignored for pattern
        SparseEntry {
            row: 0,
            col: 3,
            value: 1.0,
        },
        SparseEntry {
            row: 1,
            col: 1,
            value: 1.0,
        },
        SparseEntry {
            row: 1,
            col: 4,
            value: 1.0,
        },
        SparseEntry {
            row: 2,
            col: 2,
            value: 1.0,
        },
        SparseEntry {
            row: 3,
            col: 0,
            value: 1.0,
        },
        SparseEntry {
            row: 3,
            col: 3,
            value: 1.0,
        },
        SparseEntry {
            row: 4,
            col: 1,
            value: 1.0,
        },
        SparseEntry {
            row: 4,
            col: 4,
            value: 1.0,
        },
    ];

    let pattern_matrix = MMSparseMatrix {
        header: pattern_header,
        rows: 5,
        cols: 5,
        nnz: pattern_entries.len(),
        entries: pattern_entries,
    };

    write_sparse_matrix("pattern_matrix.mtx", &pattern_matrix)?;
    println!(
        "    Written pattern matrix (5x5) with {} non-zero positions",
        pattern_matrix.nnz
    );

    // Create an integer matrix
    println!("  Creating integer matrix...");
    let int_header = MMHeader {
        object: "matrix".to_string(),
        format: MMFormat::Coordinate,
        data_type: MMDataType::Integer,
        symmetry: MMSymmetry::General,
        comments: vec!["Integer matrix example".to_string()],
    };

    let int_entries = vec![
        SparseEntry {
            row: 0,
            col: 0,
            value: 10.0,
        },
        SparseEntry {
            row: 1,
            col: 1,
            value: 20.0,
        },
        SparseEntry {
            row: 2,
            col: 2,
            value: 30.0,
        },
        SparseEntry {
            row: 0,
            col: 2,
            value: -5.0,
        },
        SparseEntry {
            row: 2,
            col: 0,
            value: 15.0,
        },
    ];

    let int_matrix = MMSparseMatrix {
        header: int_header,
        rows: 3,
        cols: 3,
        nnz: int_entries.len(),
        entries: int_entries,
    };

    write_sparse_matrix("integer_matrix.mtx", &int_matrix)?;
    println!(
        "    Written integer matrix (3x3) with {} entries",
        int_matrix.nnz
    );

    Ok(())
}

#[allow(dead_code)]
fn demonstrate_coordinate_conversions() -> Result<()> {
    println!("\n5. Demonstrating coordinate format conversions...");

    // Create some test data in coordinate format
    let rows = Array1::from(vec![0, 1, 2, 0, 2]);
    let cols = Array1::from(vec![0, 1, 2, 2, 0]);
    let values = Array1::from(vec![1.0, 2.0, 3.0, 0.5, -1.5]);

    println!("  Original coordinate arrays:");
    println!(
        "    Rows:   {:?}",
        rows.as_slice().expect("Operation failed")
    );
    println!(
        "    Cols:   {:?}",
        cols.as_slice().expect("Operation failed")
    );
    println!(
        "    Values: {:?}",
        values.as_slice().expect("Operation failed")
    );

    // Convert to Matrix Market sparse matrix
    let header = MMHeader {
        object: "matrix".to_string(),
        format: MMFormat::Coordinate,
        data_type: MMDataType::Real,
        symmetry: MMSymmetry::General,
        comments: vec!["Converted from coordinate arrays".to_string()],
    };

    let mm_matrix = coo_to_sparse(&rows, &cols, &values, (3, 3), header);

    println!("  Converted to Matrix Market format:");
    println!("    Size: {}x{}", mm_matrix.rows, mm_matrix.cols);
    println!("    Non-zeros: {}", mm_matrix.nnz);

    // Write to file
    write_sparse_matrix("converted_matrix.mtx", &mm_matrix)?;

    // Read it back and convert back to coordinate format
    let read_matrix = read_sparse_matrix("converted_matrix.mtx")?;
    let (back_rows, back_cols, back_values) = sparse_to_coo(&read_matrix);

    println!("  Converted back to coordinate arrays:");
    println!(
        "    Rows:   {:?}",
        back_rows.as_slice().expect("Operation failed")
    );
    println!(
        "    Cols:   {:?}",
        back_cols.as_slice().expect("Operation failed")
    );
    println!(
        "    Values: {:?}",
        back_values.as_slice().expect("Operation failed")
    );

    // Verify the round-trip conversion
    let rows_match = rows == back_rows;
    let cols_match = cols == back_cols;
    let values_match = values
        .iter()
        .zip(back_values.iter())
        .all(|(&a, &b)| (a - b).abs() < 1e-10);

    println!("  Round-trip conversion verification:");
    println!("    Rows match: {}", rows_match);
    println!("    Cols match: {}", cols_match);
    println!("    Values match: {}", values_match);

    if rows_match && cols_match && values_match {
        println!("    ✓ Round-trip conversion successful!");
    } else {
        println!("    ✗ Round-trip conversion failed!");
    }

    Ok(())
}

/// Clean up example files
#[allow(dead_code)]
fn cleanup_files() {
    let files = [
        "example_sparse.mtx",
        "example_dense.mtx",
        "symmetric_matrix.mtx",
        "pattern_matrix.mtx",
        "integer_matrix.mtx",
        "converted_matrix.mtx",
    ];

    for file in &files {
        let _ = fs::remove_file(file);
    }
}