Skip to main content

math_audio_bem/core/io/
native.rs

1//! Native Rust JSON/TOML format for BEM solver configuration
2//!
3//! This module provides a clean, idiomatic Rust configuration format
4//! that supports both JSON and TOML serialization via serde.
5//!
6//! ## Example JSON Configuration
7//!
8//! ```json
9//! {
10//!     "physics": {
11//!         "frequency": 1000.0,
12//!         "speed_of_sound": 343.0,
13//!         "density": 1.21
14//!     },
15//!     "mesh": {
16//!         "nodes_file": "nodes.json",
17//!         "elements_file": "elements.json"
18//!     },
19//!     "solver": {
20//!         "method": "BiCGSTAB",
21//!         "tolerance": 1e-6,
22//!         "max_iterations": 1000
23//!     }
24//! }
25//! ```
26
27use std::fs;
28use std::path::{Path, PathBuf};
29
30use ndarray::{Array1, Array2};
31use num_complex::Complex64;
32use serde::{Deserialize, Serialize};
33
34use crate::core::types::{
35    BemMethod, BoundaryCondition, Element, ElementProperty, ElementType, PhysicsParams,
36    SolverMethod,
37};
38
39/// Native BEM configuration
40#[derive(Debug, Clone, Serialize, Deserialize)]
41pub struct BemConfig {
42    /// Problem description
43    #[serde(default)]
44    pub description: String,
45
46    /// Physics parameters
47    pub physics: PhysicsConfig,
48
49    /// Mesh configuration
50    pub mesh: MeshConfig,
51
52    /// Solver configuration
53    #[serde(default)]
54    pub solver: SolverConfig,
55
56    /// BEM method configuration
57    #[serde(default)]
58    pub bem: BemMethodConfig,
59
60    /// Boundary conditions
61    #[serde(default)]
62    pub boundary_conditions: Vec<BoundaryConditionConfig>,
63
64    /// Excitation sources
65    #[serde(default)]
66    pub sources: SourceConfig,
67
68    /// Output configuration
69    #[serde(default)]
70    pub output: OutputConfig,
71}
72
73/// Physics parameters
74#[derive(Debug, Clone, Serialize, Deserialize)]
75pub struct PhysicsConfig {
76    /// Frequency in Hz
77    pub frequency: f64,
78
79    /// Speed of sound in m/s
80    #[serde(default = "default_speed_of_sound")]
81    pub speed_of_sound: f64,
82
83    /// Medium density in kg/m³
84    #[serde(default = "default_density")]
85    pub density: f64,
86
87    /// Reference pressure in Pa
88    #[serde(default = "default_reference_pressure")]
89    pub reference_pressure: f64,
90
91    /// External problem (true) or internal (false)
92    #[serde(default = "default_true")]
93    pub external_problem: bool,
94}
95
96fn default_speed_of_sound() -> f64 {
97    343.0
98}
99fn default_density() -> f64 {
100    1.21
101}
102fn default_reference_pressure() -> f64 {
103    1.0
104}
105fn default_true() -> bool {
106    true
107}
108
109/// Mesh configuration
110#[derive(Debug, Clone, Serialize, Deserialize)]
111pub struct MeshConfig {
112    /// Path to nodes file (JSON or CSV)
113    #[serde(default)]
114    pub nodes_file: Option<PathBuf>,
115
116    /// Inline nodes data [[x, y, z], ...]
117    #[serde(default)]
118    pub nodes: Option<Vec<[f64; 3]>>,
119
120    /// Path to elements file (JSON or CSV)
121    #[serde(default)]
122    pub elements_file: Option<PathBuf>,
123
124    /// Inline elements data [[n1, n2, n3, ...], ...]
125    #[serde(default)]
126    pub elements: Option<Vec<Vec<usize>>>,
127
128    /// Symmetry planes (x, y, z)
129    #[serde(default)]
130    pub symmetry: Option<[bool; 3]>,
131
132    /// Symmetry origin
133    #[serde(default)]
134    pub symmetry_origin: Option<[f64; 3]>,
135}
136
137/// Solver configuration
138#[derive(Debug, Clone, Serialize, Deserialize, Default)]
139pub struct SolverConfig {
140    /// Iterative solver method
141    #[serde(default)]
142    pub method: SolverMethodConfig,
143
144    /// Convergence tolerance
145    #[serde(default = "default_tolerance")]
146    pub tolerance: f64,
147
148    /// Maximum iterations
149    #[serde(default = "default_max_iterations")]
150    pub max_iterations: usize,
151
152    /// Preconditioner type
153    #[serde(default)]
154    pub preconditioner: PreconditionerConfig,
155
156    /// Print progress interval (0 = no output)
157    #[serde(default)]
158    pub print_interval: usize,
159}
160
161fn default_tolerance() -> f64 {
162    1e-6
163}
164fn default_max_iterations() -> usize {
165    1000
166}
167
168/// Iterative solver method
169#[derive(Debug, Clone, Copy, Serialize, Deserialize, Default)]
170#[serde(rename_all = "lowercase")]
171pub enum SolverMethodConfig {
172    /// Direct LU factorization
173    Direct,
174    /// Conjugate Gradient Squared
175    #[default]
176    Cgs,
177    /// Bi-Conjugate Gradient Stabilized
178    #[serde(alias = "bicgstab")]
179    BiCgstab,
180}
181
182/// Preconditioner type
183#[derive(Debug, Clone, Copy, Serialize, Deserialize, Default)]
184#[serde(rename_all = "lowercase")]
185pub enum PreconditionerConfig {
186    /// No preconditioning
187    #[default]
188    None,
189    /// Diagonal (Jacobi) preconditioner
190    Diagonal,
191    /// Row scaling preconditioner
192    RowScaling,
193    /// Block diagonal preconditioner
194    BlockDiagonal,
195}
196
197/// BEM method configuration
198#[derive(Debug, Clone, Serialize, Deserialize, Default)]
199pub struct BemMethodConfig {
200    /// BEM assembly method
201    #[serde(default)]
202    pub method: BemMethodType,
203
204    /// Use Burton-Miller formulation
205    #[serde(default)]
206    pub burton_miller: bool,
207
208    /// Burton-Miller coupling parameter
209    #[serde(default = "default_coupling")]
210    pub coupling_parameter: f64,
211
212    /// Cluster size for FMM
213    #[serde(default = "default_cluster_size")]
214    pub cluster_size: usize,
215
216    /// Expansion terms for FMM
217    #[serde(default = "default_expansion_terms")]
218    pub expansion_terms: usize,
219}
220
221fn default_coupling() -> f64 {
222    1.0
223}
224fn default_cluster_size() -> usize {
225    20
226}
227fn default_expansion_terms() -> usize {
228    4
229}
230
231/// BEM assembly method type
232#[derive(Debug, Clone, Copy, Serialize, Deserialize, Default)]
233#[serde(rename_all = "lowercase")]
234pub enum BemMethodType {
235    /// Traditional BEM (dense matrices)
236    #[default]
237    Traditional,
238    /// Single-Level Fast Multipole Method
239    #[serde(alias = "slfmm")]
240    SingleLevelFmm,
241    /// Multi-Level Fast Multipole Method
242    #[serde(alias = "mlfmm")]
243    MultiLevelFmm,
244}
245
246/// Boundary condition configuration
247#[derive(Debug, Clone, Serialize, Deserialize)]
248pub struct BoundaryConditionConfig {
249    /// Element range (start, end) - inclusive
250    pub elements: (usize, usize),
251
252    /// Boundary condition type
253    #[serde(rename = "type")]
254    pub bc_type: BoundaryConditionType,
255
256    /// Value (real part)
257    pub value: f64,
258
259    /// Value (imaginary part)
260    #[serde(default)]
261    pub value_imag: f64,
262}
263
264/// Boundary condition type
265#[derive(Debug, Clone, Copy, Serialize, Deserialize)]
266#[serde(rename_all = "lowercase")]
267pub enum BoundaryConditionType {
268    /// Velocity boundary condition
269    Velocity,
270    /// Pressure boundary condition
271    Pressure,
272    /// Admittance boundary condition
273    Admittance,
274    /// Impedance boundary condition
275    Impedance,
276    /// Transfer admittance
277    TransferAdmittance,
278}
279
280/// Source configuration
281#[derive(Debug, Clone, Serialize, Deserialize, Default)]
282pub struct SourceConfig {
283    /// Plane wave sources
284    #[serde(default)]
285    pub plane_waves: Vec<PlaneWaveConfig>,
286
287    /// Point sources
288    #[serde(default)]
289    pub point_sources: Vec<PointSourceConfig>,
290}
291
292/// Plane wave source configuration
293#[derive(Debug, Clone, Serialize, Deserialize)]
294pub struct PlaneWaveConfig {
295    /// Direction vector [x, y, z] (will be normalized)
296    pub direction: [f64; 3],
297
298    /// Amplitude (complex)
299    pub amplitude: f64,
300
301    /// Phase in radians
302    #[serde(default)]
303    pub phase: f64,
304}
305
306/// Point source configuration
307#[derive(Debug, Clone, Serialize, Deserialize)]
308pub struct PointSourceConfig {
309    /// Position [x, y, z]
310    pub position: [f64; 3],
311
312    /// Amplitude (complex)
313    pub amplitude: f64,
314
315    /// Phase in radians
316    #[serde(default)]
317    pub phase: f64,
318}
319
320/// Output configuration
321#[derive(Debug, Clone, Serialize, Deserialize, Default)]
322pub struct OutputConfig {
323    /// Output directory
324    #[serde(default)]
325    pub directory: Option<PathBuf>,
326
327    /// Output pressure at evaluation points
328    #[serde(default)]
329    pub pressure: bool,
330
331    /// Output velocity at evaluation points
332    #[serde(default)]
333    pub velocity: bool,
334
335    /// Evaluation points file
336    #[serde(default)]
337    pub evaluation_points_file: Option<PathBuf>,
338
339    /// Inline evaluation points
340    #[serde(default)]
341    pub evaluation_points: Option<Vec<[f64; 3]>>,
342}
343
344/// Configuration file format
345#[derive(Debug, Clone, Copy)]
346pub enum ConfigFormat {
347    /// JSON format
348    Json,
349    /// TOML format
350    Toml,
351}
352
353impl ConfigFormat {
354    /// Detect format from file extension
355    pub fn from_path<P: AsRef<Path>>(path: P) -> Option<Self> {
356        let ext = path.as_ref().extension()?.to_str()?;
357        match ext.to_lowercase().as_str() {
358            "json" => Some(ConfigFormat::Json),
359            "toml" => Some(ConfigFormat::Toml),
360            _ => None,
361        }
362    }
363}
364
365/// Load BEM configuration from a file
366///
367/// Format is auto-detected from file extension (.json or .toml)
368pub fn load_config<P: AsRef<Path>>(path: P) -> Result<BemConfig, ConfigError> {
369    let path = path.as_ref();
370    let content = fs::read_to_string(path)?;
371
372    let format = ConfigFormat::from_path(path)
373        .ok_or_else(|| ConfigError::UnsupportedFormat(path.display().to_string()))?;
374
375    parse_config(&content, format)
376}
377
378/// Parse BEM configuration from a string
379pub fn parse_config(content: &str, format: ConfigFormat) -> Result<BemConfig, ConfigError> {
380    match format {
381        ConfigFormat::Json => {
382            serde_json::from_str(content).map_err(|e| ConfigError::ParseError(e.to_string()))
383        }
384        ConfigFormat::Toml => {
385            toml::from_str(content).map_err(|e| ConfigError::ParseError(e.to_string()))
386        }
387    }
388}
389
390/// Save BEM configuration to a file
391pub fn save_config<P: AsRef<Path>>(config: &BemConfig, path: P) -> Result<(), ConfigError> {
392    let path = path.as_ref();
393    let format = ConfigFormat::from_path(path)
394        .ok_or_else(|| ConfigError::UnsupportedFormat(path.display().to_string()))?;
395
396    let content = serialize_config(config, format)?;
397    fs::write(path, content)?;
398    Ok(())
399}
400
401/// Serialize BEM configuration to a string
402pub fn serialize_config(config: &BemConfig, format: ConfigFormat) -> Result<String, ConfigError> {
403    match format {
404        ConfigFormat::Json => serde_json::to_string_pretty(config)
405            .map_err(|e| ConfigError::SerializeError(e.to_string())),
406        ConfigFormat::Toml => {
407            toml::to_string_pretty(config).map_err(|e| ConfigError::SerializeError(e.to_string()))
408        }
409    }
410}
411
412/// Configuration error types
413#[derive(Debug, thiserror::Error)]
414pub enum ConfigError {
415    /// IO error
416    #[error("IO error: {0}")]
417    Io(#[from] std::io::Error),
418
419    /// Parse error
420    #[error("Parse error: {0}")]
421    ParseError(String),
422
423    /// Serialize error
424    #[error("Serialize error: {0}")]
425    SerializeError(String),
426
427    /// Unsupported format
428    #[error("Unsupported format: {0}")]
429    UnsupportedFormat(String),
430
431    /// Missing required field
432    #[error("Missing required field: {0}")]
433    MissingField(String),
434}
435
436impl BemConfig {
437    /// Create PhysicsParams from configuration
438    pub fn to_physics_params(&self) -> PhysicsParams {
439        PhysicsParams::new(
440            self.physics.frequency,
441            self.physics.speed_of_sound,
442            self.physics.density,
443            self.physics.external_problem,
444        )
445    }
446
447    /// Get BEM method from configuration
448    pub fn bem_method(&self) -> BemMethod {
449        match self.bem.method {
450            BemMethodType::Traditional => BemMethod::Traditional,
451            BemMethodType::SingleLevelFmm => BemMethod::SingleLevelFmm,
452            BemMethodType::MultiLevelFmm => BemMethod::MultiLevelFmm,
453        }
454    }
455
456    /// Get solver method from configuration
457    pub fn solver_method(&self) -> SolverMethod {
458        match self.solver.method {
459            SolverMethodConfig::Direct => SolverMethod::Direct,
460            SolverMethodConfig::Cgs => SolverMethod::Cgs,
461            SolverMethodConfig::BiCgstab => SolverMethod::BiCgstab,
462        }
463    }
464
465    /// Load nodes from configuration
466    pub fn load_nodes(&self, base_dir: &Path) -> Result<Array2<f64>, ConfigError> {
467        if let Some(ref nodes) = self.mesh.nodes {
468            // Inline nodes
469            let n = nodes.len();
470            let mut arr = Array2::zeros((n, 3));
471            for (i, node) in nodes.iter().enumerate() {
472                arr[[i, 0]] = node[0];
473                arr[[i, 1]] = node[1];
474                arr[[i, 2]] = node[2];
475            }
476            return Ok(arr);
477        }
478
479        if let Some(ref file) = self.mesh.nodes_file {
480            let path = base_dir.join(file);
481            let content = fs::read_to_string(&path)?;
482
483            // Try JSON first
484            if path.extension().is_some_and(|e| e == "json") {
485                let nodes: Vec<[f64; 3]> = serde_json::from_str(&content)
486                    .map_err(|e| ConfigError::ParseError(e.to_string()))?;
487                let n = nodes.len();
488                let mut arr = Array2::zeros((n, 3));
489                for (i, node) in nodes.iter().enumerate() {
490                    arr[[i, 0]] = node[0];
491                    arr[[i, 1]] = node[1];
492                    arr[[i, 2]] = node[2];
493                }
494                return Ok(arr);
495            }
496
497            // Try CSV
498            let nodes = parse_csv_nodes(&content)?;
499            return Ok(nodes);
500        }
501
502        Err(ConfigError::MissingField(
503            "mesh.nodes or mesh.nodes_file".to_string(),
504        ))
505    }
506
507    /// Load elements from configuration
508    pub fn load_elements(&self, base_dir: &Path) -> Result<Vec<Element>, ConfigError> {
509        if let Some(ref elements) = self.mesh.elements {
510            return Ok(elements_from_connectivity(elements));
511        }
512
513        if let Some(ref file) = self.mesh.elements_file {
514            let path = base_dir.join(file);
515            let content = fs::read_to_string(&path)?;
516
517            // Try JSON first
518            if path.extension().is_some_and(|e| e == "json") {
519                let connectivity: Vec<Vec<usize>> = serde_json::from_str(&content)
520                    .map_err(|e| ConfigError::ParseError(e.to_string()))?;
521                return Ok(elements_from_connectivity(&connectivity));
522            }
523
524            // Try CSV
525            let elements = parse_csv_elements(&content)?;
526            return Ok(elements);
527        }
528
529        Err(ConfigError::MissingField(
530            "mesh.elements or mesh.elements_file".to_string(),
531        ))
532    }
533}
534
535/// Parse CSV nodes (x y z per line)
536fn parse_csv_nodes(content: &str) -> Result<Array2<f64>, ConfigError> {
537    let mut nodes = Vec::new();
538
539    for line in content.lines() {
540        let line = line.trim();
541        if line.is_empty() || line.starts_with('#') {
542            continue;
543        }
544
545        let values: Vec<f64> = line
546            .split([',', ' ', '\t'])
547            .filter_map(|s| s.trim().parse().ok())
548            .collect();
549
550        if values.len() >= 3 {
551            nodes.push([values[0], values[1], values[2]]);
552        }
553    }
554
555    let n = nodes.len();
556    let mut arr = Array2::zeros((n, 3));
557    for (i, node) in nodes.iter().enumerate() {
558        arr[[i, 0]] = node[0];
559        arr[[i, 1]] = node[1];
560        arr[[i, 2]] = node[2];
561    }
562
563    Ok(arr)
564}
565
566/// Parse CSV elements (connectivity per line)
567fn parse_csv_elements(content: &str) -> Result<Vec<Element>, ConfigError> {
568    let mut connectivity = Vec::new();
569
570    for line in content.lines() {
571        let line = line.trim();
572        if line.is_empty() || line.starts_with('#') {
573            continue;
574        }
575
576        let values: Vec<usize> = line
577            .split([',', ' ', '\t'])
578            .filter_map(|s| s.trim().parse().ok())
579            .collect();
580
581        if values.len() >= 3 {
582            connectivity.push(values);
583        }
584    }
585
586    Ok(elements_from_connectivity(&connectivity))
587}
588
589/// Create elements from connectivity data
590fn elements_from_connectivity(connectivity: &[Vec<usize>]) -> Vec<Element> {
591    connectivity
592        .iter()
593        .enumerate()
594        .map(|(idx, conn)| {
595            let element_type = if conn.len() == 3 {
596                ElementType::Tri3
597            } else {
598                ElementType::Quad4
599            };
600
601            Element {
602                connectivity: conn.clone(),
603                element_type,
604                property: ElementProperty::Surface,
605                normal: Array1::zeros(3),
606                node_normals: Array2::zeros((element_type.num_nodes(), 3)),
607                center: Array1::zeros(3),
608                area: 0.0,
609                boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(0.0, 0.0)]),
610                group: 0,
611                dof_addresses: vec![idx],
612            }
613        })
614        .collect()
615}
616
617#[cfg(test)]
618mod tests {
619    use super::*;
620
621    const SAMPLE_JSON: &str = r#"{
622        "description": "Test BEM problem",
623        "physics": {
624            "frequency": 1000.0,
625            "speed_of_sound": 343.0,
626            "density": 1.21
627        },
628        "mesh": {
629            "nodes": [[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [0.0, 1.0, 0.0]],
630            "elements": [[0, 1, 2]]
631        },
632        "solver": {
633            "method": "bicgstab",
634            "tolerance": 1e-8,
635            "max_iterations": 500
636        },
637        "bem": {
638            "method": "traditional",
639            "burton_miller": true
640        },
641        "boundary_conditions": [
642            {
643                "elements": [0, 0],
644                "type": "velocity",
645                "value": 1.0
646            }
647        ],
648        "sources": {
649            "plane_waves": [
650                {
651                    "direction": [0.0, 0.0, -1.0],
652                    "amplitude": 1.0,
653                    "phase": 0.0
654                }
655            ]
656        }
657    }"#;
658
659    const SAMPLE_TOML: &str = r#"
660description = "Test BEM problem"
661
662[physics]
663frequency = 1000.0
664speed_of_sound = 343.0
665density = 1.21
666
667[mesh]
668nodes = [[0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [0.0, 1.0, 0.0]]
669elements = [[0, 1, 2]]
670
671[solver]
672method = "bicgstab"
673tolerance = 1e-8
674max_iterations = 500
675
676[bem]
677method = "traditional"
678burton_miller = true
679
680[[boundary_conditions]]
681elements = [0, 0]
682type = "velocity"
683value = 1.0
684
685[sources]
686[[sources.plane_waves]]
687direction = [0.0, 0.0, -1.0]
688amplitude = 1.0
689phase = 0.0
690"#;
691
692    #[test]
693    fn test_parse_json() {
694        let config = parse_config(SAMPLE_JSON, ConfigFormat::Json).unwrap();
695
696        assert_eq!(config.description, "Test BEM problem");
697        assert!((config.physics.frequency - 1000.0).abs() < 0.01);
698        assert!((config.physics.speed_of_sound - 343.0).abs() < 0.01);
699        assert_eq!(config.mesh.nodes.as_ref().unwrap().len(), 3);
700        assert_eq!(config.mesh.elements.as_ref().unwrap().len(), 1);
701        assert!(matches!(config.solver.method, SolverMethodConfig::BiCgstab));
702        assert!(config.bem.burton_miller);
703        assert_eq!(config.boundary_conditions.len(), 1);
704        assert_eq!(config.sources.plane_waves.len(), 1);
705    }
706
707    #[test]
708    fn test_parse_toml() {
709        let config = parse_config(SAMPLE_TOML, ConfigFormat::Toml).unwrap();
710
711        assert_eq!(config.description, "Test BEM problem");
712        assert!((config.physics.frequency - 1000.0).abs() < 0.01);
713        assert!((config.physics.speed_of_sound - 343.0).abs() < 0.01);
714        assert_eq!(config.mesh.nodes.as_ref().unwrap().len(), 3);
715    }
716
717    #[test]
718    fn test_to_physics_params() {
719        let config = parse_config(SAMPLE_JSON, ConfigFormat::Json).unwrap();
720        let physics = config.to_physics_params();
721
722        assert!((physics.frequency - 1000.0).abs() < 0.01);
723        assert!((physics.speed_of_sound - 343.0).abs() < 0.01);
724        assert!((physics.density - 1.21).abs() < 0.01);
725    }
726
727    #[test]
728    fn test_load_inline_nodes() {
729        let config = parse_config(SAMPLE_JSON, ConfigFormat::Json).unwrap();
730        let nodes = config.load_nodes(Path::new(".")).unwrap();
731
732        assert_eq!(nodes.nrows(), 3);
733        assert_eq!(nodes.ncols(), 3);
734        assert!((nodes[[0, 0]] - 0.0).abs() < 1e-10);
735        assert!((nodes[[1, 0]] - 1.0).abs() < 1e-10);
736    }
737
738    #[test]
739    fn test_load_inline_elements() {
740        let config = parse_config(SAMPLE_JSON, ConfigFormat::Json).unwrap();
741        let elements = config.load_elements(Path::new(".")).unwrap();
742
743        assert_eq!(elements.len(), 1);
744        assert!(matches!(elements[0].element_type, ElementType::Tri3));
745    }
746
747    #[test]
748    fn test_serialize_json() {
749        let config = parse_config(SAMPLE_JSON, ConfigFormat::Json).unwrap();
750        let serialized = serialize_config(&config, ConfigFormat::Json).unwrap();
751
752        // Verify we can parse it back
753        let reparsed = parse_config(&serialized, ConfigFormat::Json).unwrap();
754        assert_eq!(reparsed.description, config.description);
755    }
756
757    #[test]
758    fn test_serialize_toml() {
759        let config = parse_config(SAMPLE_JSON, ConfigFormat::Json).unwrap();
760        let serialized = serialize_config(&config, ConfigFormat::Toml).unwrap();
761
762        // Verify we can parse it back
763        let reparsed = parse_config(&serialized, ConfigFormat::Toml).unwrap();
764        assert_eq!(reparsed.description, config.description);
765    }
766
767    #[test]
768    fn test_bem_method_conversion() {
769        let config = parse_config(SAMPLE_JSON, ConfigFormat::Json).unwrap();
770        assert!(matches!(config.bem_method(), BemMethod::Traditional));
771
772        let json_fmm = r#"{
773            "physics": {"frequency": 1000.0},
774            "mesh": {"nodes": [], "elements": []},
775            "bem": {"method": "slfmm"}
776        }"#;
777        let config_fmm = parse_config(json_fmm, ConfigFormat::Json).unwrap();
778        assert!(matches!(config_fmm.bem_method(), BemMethod::SingleLevelFmm));
779    }
780
781    #[test]
782    fn test_solver_method_conversion() {
783        let config = parse_config(SAMPLE_JSON, ConfigFormat::Json).unwrap();
784        assert!(matches!(config.solver_method(), SolverMethod::BiCgstab));
785    }
786
787    #[test]
788    fn test_parse_csv_nodes() {
789        let csv = "0.0 0.0 0.0\n1.0 0.0 0.0\n0.5 1.0 0.0";
790        let nodes = parse_csv_nodes(csv).unwrap();
791
792        assert_eq!(nodes.nrows(), 3);
793        assert!((nodes[[1, 0]] - 1.0).abs() < 1e-10);
794    }
795
796    #[test]
797    fn test_parse_csv_elements() {
798        let csv = "0, 1, 2\n1, 2, 3";
799        let elements = parse_csv_elements(csv).unwrap();
800
801        assert_eq!(elements.len(), 2);
802        assert_eq!(elements[0].connectivity, vec![0, 1, 2]);
803    }
804}