Skip to main content

causal_triangulations/
lib.rs

1#![allow(clippy::multiple_crate_versions)]
2#![warn(missing_docs)]
3
4//! Causal Dynamical Triangulations library for quantum gravity simulations.
5//!
6//! This library implements Causal Dynamical Triangulations (CDT) in 2D, providing
7//! the necessary tools for Monte Carlo simulations of discrete spacetime geometries.
8//!
9//! # Key Features
10//!
11//! - Integration with delaunay crate for proper Delaunay triangulations
12//! - 2D Regge Action calculation for CDT
13//! - Standard ergodic moves: (2,2), (1,3), and edge flips
14//! - Metropolis-Hastings algorithm with configurable temperature
15//! - Comprehensive statistics and measurement collection
16//!
17//! # Example
18//!
19//! ```rust,no_run
20//! // Example would require command line arguments, so we skip execution
21//! use causal_triangulations::{CdtConfig, run_simulation};
22//! // CdtConfig requires configuration, so this is marked no_run
23//! ```
24
25// Module declarations (avoiding mod.rs files)
26/// Configuration management for CDT simulations.
27pub mod config;
28
29/// Error types for the CDT library.
30pub mod errors;
31
32/// Utility functions for random number generation and mathematical operations.
33pub mod util;
34
35/// Geometry abstraction layer for CDT simulations.
36///
37/// This module provides trait-based geometry operations that isolate CDT algorithms
38/// from specific geometry implementations.
39pub mod geometry {
40    /// CDT-agnostic mesh data structures.
41    pub mod mesh;
42    /// High-level triangulation operations.
43    pub mod operations;
44    /// Core geometry traits for CDT abstraction.
45    pub mod traits;
46
47    /// Geometry backend implementations.
48    pub mod backends {
49        /// Delaunay backend - wraps the delaunay crate.
50        pub mod delaunay;
51
52        /// Mock backend for testing.
53        pub mod mock;
54    }
55
56    // Type aliases for common backend combinations
57    /// 2D Delaunay backend (most common configuration).
58    ///
59    /// Uses `f64` coordinates with `()` vertex data and `i32` cell data.
60    /// CDT metadata is tracked separately by [`CdtTriangulation`](crate::cdt::triangulation::CdtTriangulation).
61    pub type DelaunayBackend2D = backends::delaunay::DelaunayBackend<(), i32, 2>;
62
63    /// Default backend type for 2D CDT simulations
64    pub type DefaultBackend = DelaunayBackend2D;
65
66    /// Convenient alias for CDT triangulations using the default backend
67    pub type CdtTriangulation2D = crate::cdt::triangulation::CdtTriangulation<DefaultBackend>;
68}
69
70/// Causal Dynamical Triangulations implementation modules.
71pub mod cdt {
72    /// Action calculation for CDT simulations.
73    pub mod action;
74    /// Ergodic moves for triangulation modifications.
75    pub mod ergodic_moves;
76    /// Metropolis-Hastings algorithm implementation.
77    pub mod metropolis;
78    /// CDT triangulation wrapper.
79    pub mod triangulation;
80}
81
82// Re-exports for convenience
83pub use cdt::action::{ActionConfig, compute_regge_action};
84pub use cdt::ergodic_moves::{ErgodicsSystem, MoveResult, MoveType};
85pub use cdt::metropolis::{MetropolisAlgorithm, MetropolisConfig, SimulationResultsBackend};
86pub use config::{CdtConfig, TestConfig};
87pub use errors::{CdtError, CdtResult};
88
89use cdt::metropolis::Measurement;
90use std::time::Duration;
91
92// Trait-based triangulation (recommended)
93pub use cdt::triangulation::CdtTriangulation;
94
95/// Runs a CDT simulation with the specified configuration.
96///
97/// This function uses the trait-based geometry backend system, which provides
98/// better abstraction and testability compared to legacy approaches.
99///
100/// # Arguments
101///
102/// * `config` - Configuration parameters for the triangulation/simulation
103///
104/// # Returns
105///
106/// A `SimulationResults` struct containing the results of the simulation.
107///
108/// # Errors
109///
110/// Returns [`CdtError::InvalidParameters`] if the configuration fails validation
111/// (e.g., invalid measurement frequency, inconsistent parameters).
112/// Returns [`CdtError::UnsupportedDimension`] if an unsupported dimension (not 2D) is specified.
113/// Returns triangulation generation errors from the underlying triangulation creation.
114pub fn run_simulation(config: &CdtConfig) -> CdtResult<cdt::metropolis::SimulationResultsBackend> {
115    // Validate configuration early to fail fast with clear error messages
116    if let Err(err) = config.validate() {
117        return Err(CdtError::InvalidParameters(err));
118    }
119
120    let vertices = config.vertices;
121    let timeslices = config.timeslices;
122
123    if let Some(dim) = config.dimension
124        && dim != 2
125    {
126        return Err(CdtError::UnsupportedDimension(dim.into()));
127    }
128
129    log::info!("Dimensionality: {}", config.dimension.unwrap_or(2));
130    log::info!("Number of vertices: {vertices}");
131    log::info!("Number of timeslices: {timeslices}");
132    log::info!("Using trait-based backend system");
133
134    // Create initial triangulation
135    let triangulation = CdtTriangulation::from_random_points(vertices, timeslices, 2)?;
136
137    log::info!(
138        "Triangulation created with {} vertices, {} edges, {} faces",
139        triangulation.vertex_count(),
140        triangulation.edge_count(),
141        triangulation.face_count()
142    );
143
144    if config.simulate {
145        // Run full CDT simulation with new backend
146        let metropolis_config = config.to_metropolis_config();
147        let action_config = config.to_action_config();
148
149        let mut algorithm = MetropolisAlgorithm::new(metropolis_config, action_config);
150        let results = algorithm.run(triangulation);
151
152        log::info!("Simulation Results:");
153        log::info!(
154            "  Acceptance rate: {:.2}%",
155            results.acceptance_rate() * 100.0
156        );
157        log::info!("  Average action: {:.3}", results.average_action());
158
159        Ok(results)
160    } else {
161        // Just return basic simulation results with the triangulation
162        let initial_action = config.to_action_config().calculate_action(
163            u32::try_from(triangulation.vertex_count()).unwrap_or_default(),
164            u32::try_from(triangulation.edge_count()).unwrap_or_default(),
165            u32::try_from(triangulation.face_count()).unwrap_or_default(),
166        );
167
168        Ok(cdt::metropolis::SimulationResultsBackend {
169            config: config.to_metropolis_config(),
170            action_config: config.to_action_config(),
171            steps: vec![],
172            measurements: vec![Measurement {
173                step: 0,
174                action: initial_action,
175                vertices: u32::try_from(triangulation.vertex_count()).unwrap_or_default(),
176                edges: u32::try_from(triangulation.edge_count()).unwrap_or_default(),
177                triangles: u32::try_from(triangulation.face_count()).unwrap_or_default(),
178            }],
179            elapsed_time: Duration::from_millis(0),
180            triangulation,
181        })
182    }
183}
184
185#[cfg(test)]
186mod lib_tests {
187    use super::*;
188    use approx::assert_relative_eq;
189
190    fn create_test_config() -> CdtConfig {
191        CdtConfig {
192            dimension: Some(2),
193            vertices: 32,
194            timeslices: 3,
195            temperature: 1.0,
196            steps: 10,
197            thermalization_steps: 5,
198            measurement_frequency: 2,
199            coupling_0: 1.0,
200            coupling_2: 1.0,
201            cosmological_constant: 0.1,
202            simulate: false,
203        }
204    }
205
206    #[test]
207    fn test_run_simulation() {
208        let config = create_test_config();
209        assert!(config.dimension.is_some());
210        let results = run_simulation(&config).expect("Failed to run triangulation");
211        assert!(results.triangulation.face_count() > 0);
212        assert!(!results.measurements.is_empty());
213    }
214
215    #[test]
216    fn triangulation_contains_triangles() {
217        let config = create_test_config();
218        let results = run_simulation(&config).expect("Failed to run triangulation");
219        // Check that we have some triangles
220        assert!(results.triangulation.face_count() > 0);
221    }
222
223    #[test]
224    fn test_config_validation_invalid_measurement_frequency() {
225        let mut config = create_test_config();
226        config.measurement_frequency = 0;
227
228        let result = run_simulation(&config);
229        assert!(result.is_err(), "Should reject zero measurement frequency");
230
231        if let Err(CdtError::InvalidParameters(msg)) = result {
232            assert!(
233                msg.contains("Measurement frequency must be positive"),
234                "Error message should mention measurement frequency: {msg}"
235            );
236        } else {
237            panic!("Expected InvalidParameters error");
238        }
239    }
240
241    #[test]
242    fn test_config_validation_measurement_frequency_too_large() {
243        let mut config = create_test_config();
244        config.steps = 100;
245        config.measurement_frequency = 200; // Greater than steps
246
247        let result = run_simulation(&config);
248        assert!(
249            result.is_err(),
250            "Should reject measurement frequency greater than steps"
251        );
252
253        if let Err(CdtError::InvalidParameters(msg)) = result {
254            assert!(
255                msg.contains("cannot be greater than total steps"),
256                "Error message should mention steps constraint: {msg}"
257            );
258        } else {
259            panic!("Expected InvalidParameters error");
260        }
261    }
262
263    #[test]
264    fn test_config_validation_invalid_vertices() {
265        let mut config = create_test_config();
266        config.vertices = 2; // Less than minimum of 3
267
268        let result = run_simulation(&config);
269        assert!(result.is_err(), "Should reject too few vertices");
270
271        if let Err(CdtError::InvalidParameters(msg)) = result {
272            assert!(
273                msg.contains("vertices must be at least 3"),
274                "Error message should mention vertex constraint: {msg}"
275            );
276        } else {
277            panic!("Expected InvalidParameters error");
278        }
279    }
280
281    #[test]
282    fn test_config_validation_negative_temperature() {
283        let mut config = create_test_config();
284        config.temperature = -1.0;
285
286        let result = run_simulation(&config);
287        assert!(result.is_err(), "Should reject negative temperature");
288
289        if let Err(CdtError::InvalidParameters(msg)) = result {
290            assert!(
291                msg.contains("Temperature must be positive"),
292                "Error message should mention temperature constraint: {msg}"
293            );
294        } else {
295            panic!("Expected InvalidParameters error");
296        }
297    }
298
299    #[test]
300    fn test_config_conversions() {
301        let config = create_test_config();
302
303        let metropolis_config = config.to_metropolis_config();
304        assert_relative_eq!(metropolis_config.temperature, 1.0);
305        assert_eq!(metropolis_config.steps, 10);
306
307        let action_config = config.to_action_config();
308        assert_relative_eq!(action_config.coupling_0, 1.0);
309        assert_relative_eq!(action_config.coupling_2, 1.0);
310        assert_relative_eq!(action_config.cosmological_constant, 0.1);
311    }
312}