Skip to main content

caustic/tooling/core/init/
cosmological.rs

1//! Zel'dovich approximation initial conditions for cosmological structure formation.
2//! f is a 3D sheet in 6D phase space (cold dark matter).
3
4use super::super::types::PhaseSpaceSnapshot;
5use super::domain::Domain;
6use rust_decimal::prelude::ToPrimitive;
7
8/// Cosmological parameters for Friedmann background.
9pub struct CosmologyParams {
10    pub h0: f64,
11    pub omega_m: f64,
12    pub omega_lambda: f64,
13    pub a_init: f64,
14}
15
16/// 1D matter power spectrum P(k) for seeding perturbations.
17pub struct PowerSpectrum {
18    /// (k, P(k)) pairs.
19    pub values: Vec<(f64, f64)>,
20}
21
22/// Zel'dovich pancake IC: cold dark matter sheet.
23/// f(x,v,t₀) = ρ̄·δ³(v − v₀(x)) where v₀ is the Zel'dovich velocity field.
24pub struct ZeldovichIC {
25    pub mean_density: f64,
26    pub h0: f64,
27    pub omega_m: f64,
28    pub omega_lambda: f64,
29    pub scale_factor_init: f64,
30    pub random_seed: u64,
31}
32
33impl ZeldovichIC {
34    pub fn new(mean_density: f64, cosmology: CosmologyParams, seed: u64) -> Self {
35        Self {
36            mean_density,
37            h0: cosmology.h0,
38            omega_m: cosmology.omega_m,
39            omega_lambda: cosmology.omega_lambda,
40            scale_factor_init: cosmology.a_init,
41            random_seed: seed,
42        }
43    }
44
45    /// Create a simple single-mode Zel'dovich pancake for testing.
46    ///
47    /// v₀(x) = -A·sin(k·x₁)
48    /// f(x,v) = (ρ̄/S) · exp(-(v - v₀(x))² / (2σ_v²))
49    ///
50    /// where S normalises the Gaussian to unity over velocity space.
51    pub fn new_single_mode(
52        mean_density: f64,
53        amplitude: f64,
54        wavenumber: f64,
55        sigma_v: f64,
56    ) -> ZeldovichSingleMode {
57        ZeldovichSingleMode {
58            mean_density,
59            amplitude,
60            wavenumber,
61            sigma_v,
62        }
63    }
64
65    /// Generate the Zel'dovich displacement field s(q) from P(k) via FFT.
66    ///
67    /// Uses `Xorshift64` RNG with Box-Muller transform for Gaussian random
68    /// field generation in Fourier space, then IFFT to real space.
69    pub fn displacement_field(&self, domain: &Domain) -> [Vec<f64>; 3] {
70        use crate::tooling::core::algos::aca::Xorshift64;
71        use rustfft::num_complex::Complex64;
72        use std::f64::consts::PI;
73
74        let nx = domain.spatial_res.x1 as usize;
75        let ny = domain.spatial_res.x2 as usize;
76        let nz = domain.spatial_res.x3 as usize;
77        let n_total = nx * ny * nz;
78
79        // Box lengths: domain is [-L, L], so total length = 2L
80        let lx = 2.0 * domain.spatial.x1.to_f64().unwrap();
81        let ly = 2.0 * domain.spatial.x2.to_f64().unwrap();
82        let lz = 2.0 * domain.spatial.x3.to_f64().unwrap();
83        let volume = lx * ly * lz;
84
85        let dx_grid = lx / nx as f64;
86        let dy_grid = ly / ny as f64;
87        let dz_grid = lz / nz as f64;
88
89        // Nyquist wavenumbers
90        let kn_x = PI / dx_grid;
91        let kn_y = PI / dy_grid;
92        let kn_z = PI / dz_grid;
93        let kn2_avg = (kn_x * kn_x + kn_y * kn_y + kn_z * kn_z) / 3.0;
94
95        // Mean spacing and target RMS displacement
96        let mean_spacing = (dx_grid * dy_grid * dz_grid).cbrt();
97        let target_rms = 0.1 * mean_spacing;
98
99        // RNG + Box-Muller helper
100        let mut rng = Xorshift64::new(self.random_seed);
101
102        fn box_muller(rng: &mut Xorshift64) -> (f64, f64) {
103            let u1 = (rng.next_u64() as f64 + 1.0) / (u64::MAX as f64 + 1.0);
104            let u2 = (rng.next_u64() as f64 + 1.0) / (u64::MAX as f64 + 1.0);
105            let r = (-2.0 * u1.ln()).sqrt();
106            let theta = 2.0 * PI * u2;
107            (r * theta.cos(), r * theta.sin())
108        }
109
110        // Generate δ(k) in Fourier space with Harrison-Zel'dovich-like spectrum
111        // P(k) = A * |k| * exp(-|k|² / (2 k_N²))
112        // Amplitude A is calibrated after to get target_rms displacement.
113        let mut delta_k = vec![Complex64::new(0.0, 0.0); n_total];
114
115        // First pass: generate with A=1, then rescale
116        for i0 in 0..nx {
117            let kx_idx = if i0 <= nx / 2 {
118                i0 as f64
119            } else {
120                i0 as f64 - nx as f64
121            };
122            let kx = 2.0 * PI * kx_idx / lx;
123            for i1 in 0..ny {
124                let ky_idx = if i1 <= ny / 2 {
125                    i1 as f64
126                } else {
127                    i1 as f64 - ny as f64
128                };
129                let ky = 2.0 * PI * ky_idx / ly;
130                for i2 in 0..nz {
131                    let kz_idx = if i2 <= nz / 2 {
132                        i2 as f64
133                    } else {
134                        i2 as f64 - nz as f64
135                    };
136                    let kz = 2.0 * PI * kz_idx / lz;
137
138                    let k2 = kx * kx + ky * ky + kz * kz;
139                    let idx = i0 * ny * nz + i1 * nz + i2;
140
141                    if k2 == 0.0 {
142                        // Skip DC mode
143                        delta_k[idx] = Complex64::new(0.0, 0.0);
144                        continue;
145                    }
146
147                    let k_mag = k2.sqrt();
148                    // P(k) = |k| * exp(-k² / (2 k_N²))  (A=1 initially)
149                    let pk = k_mag * (-k2 / (2.0 * kn2_avg)).exp();
150                    let amplitude = (pk / volume).sqrt();
151
152                    let (z1, z2) = box_muller(&mut rng);
153                    delta_k[idx] = Complex64::new(z1 * amplitude, z2 * amplitude);
154                }
155            }
156        }
157
158        // Enforce Hermitian symmetry: δ(-k) = conj(δ(k))
159        for i0 in 0..nx {
160            let j0 = if i0 == 0 { 0 } else { nx - i0 };
161            for i1 in 0..ny {
162                let j1 = if i1 == 0 { 0 } else { ny - i1 };
163                for i2 in 0..nz {
164                    let j2 = if i2 == 0 { 0 } else { nz - i2 };
165
166                    let idx_pos = i0 * ny * nz + i1 * nz + i2;
167                    let idx_neg = j0 * ny * nz + j1 * nz + j2;
168
169                    if idx_pos < idx_neg {
170                        // Set conjugate pair
171                        let val = delta_k[idx_pos];
172                        delta_k[idx_neg] = val.conj();
173                    } else if idx_pos == idx_neg {
174                        // Self-conjugate mode: must be real
175                        delta_k[idx_pos] = Complex64::new(delta_k[idx_pos].re, 0.0);
176                    }
177                }
178            }
179        }
180
181        // Compute displacement in k-space: s_k = i·k / |k|² · δ(k)
182        let mut sx_k = vec![Complex64::new(0.0, 0.0); n_total];
183        let mut sy_k = vec![Complex64::new(0.0, 0.0); n_total];
184        let mut sz_k = vec![Complex64::new(0.0, 0.0); n_total];
185
186        for i0 in 0..nx {
187            let kx_idx = if i0 <= nx / 2 {
188                i0 as f64
189            } else {
190                i0 as f64 - nx as f64
191            };
192            let kx = 2.0 * PI * kx_idx / lx;
193            for i1 in 0..ny {
194                let ky_idx = if i1 <= ny / 2 {
195                    i1 as f64
196                } else {
197                    i1 as f64 - ny as f64
198                };
199                let ky = 2.0 * PI * ky_idx / ly;
200                for i2 in 0..nz {
201                    let kz_idx = if i2 <= nz / 2 {
202                        i2 as f64
203                    } else {
204                        i2 as f64 - nz as f64
205                    };
206                    let kz = 2.0 * PI * kz_idx / lz;
207
208                    let k2 = kx * kx + ky * ky + kz * kz;
209                    let idx = i0 * ny * nz + i1 * nz + i2;
210
211                    if k2 == 0.0 {
212                        continue;
213                    }
214
215                    // s_k = i·k/|k|² · δ(k)
216                    let dk = delta_k[idx];
217                    let i_dk = Complex64::new(-dk.im, dk.re); // i * δ(k)
218                    sx_k[idx] = i_dk * (kx / k2);
219                    sy_k[idx] = i_dk * (ky / k2);
220                    sz_k[idx] = i_dk * (kz / k2);
221                }
222            }
223        }
224
225        // 3D IFFT each component
226        let sx = ifft_3d(&sx_k, [nx, ny, nz]);
227        let sy = ifft_3d(&sy_k, [nx, ny, nz]);
228        let sz = ifft_3d(&sz_k, [nx, ny, nz]);
229
230        // Rescale to achieve target RMS displacement
231        let rms = {
232            let sum_sq: f64 = sx
233                .iter()
234                .chain(sy.iter())
235                .chain(sz.iter())
236                .map(|v| v * v)
237                .sum();
238            (sum_sq / (3.0 * n_total as f64)).sqrt()
239        };
240
241        if rms > 0.0 {
242            let scale = target_rms / rms;
243            let sx: Vec<f64> = sx.iter().map(|v| v * scale).collect();
244            let sy: Vec<f64> = sy.iter().map(|v| v * scale).collect();
245            let sz: Vec<f64> = sz.iter().map(|v| v * scale).collect();
246            [sx, sy, sz]
247        } else {
248            [sx, sy, sz]
249        }
250    }
251
252    /// Zel'dovich velocity v₀(q) = a·H(a)·f(a)·s(q) (growing mode).
253    ///
254    /// f(a) ≈ Ω_m(a)^{0.55} is the logarithmic growth rate.
255    pub fn velocity_field(&self, domain: &Domain) -> [Vec<f64>; 3] {
256        let [sx, sy, sz] = self.displacement_field(domain);
257
258        let a = self.scale_factor_init;
259        let a3 = a * a * a;
260
261        // E(a) = H(a)/H0
262        let e2 = self.omega_m / a3 + self.omega_lambda;
263        let h_a = self.h0 * e2.sqrt(); // H(a)
264
265        // Ω_m(a) at this scale factor
266        let om_a = self.omega_m / (a3 * e2);
267
268        // Logarithmic growth rate f ≈ Ω_m(a)^{0.55}
269        let f_growth = om_a.powf(0.55);
270
271        // v(q) = a · H(a) · f(a) · s(q)
272        let factor = a * h_a * f_growth;
273
274        let vx: Vec<f64> = sx.iter().map(|s| s * factor).collect();
275        let vy: Vec<f64> = sy.iter().map(|s| s * factor).collect();
276        let vz: Vec<f64> = sz.iter().map(|s| s * factor).collect();
277
278        [vx, vy, vz]
279    }
280
281    /// Sample onto 6D grid as a thin Gaussian in velocity centred on v₀(x).
282    ///
283    /// f(x,v) = ρ̄ / ((2π)^{3/2} σ_v³) · exp(-|v - v₀(x)|² / (2σ_v²))
284    /// where σ_v = 0.3 · max(Δv) keeps the distribution cold but resolved.
285    pub fn sample_on_grid(&self, domain: &Domain) -> PhaseSpaceSnapshot {
286        use std::f64::consts::PI;
287
288        let [vx0, vy0, vz0] = self.velocity_field(domain);
289
290        let nx1 = domain.spatial_res.x1 as usize;
291        let nx2 = domain.spatial_res.x2 as usize;
292        let nx3 = domain.spatial_res.x3 as usize;
293        let nv1 = domain.velocity_res.v1 as usize;
294        let nv2 = domain.velocity_res.v2 as usize;
295        let nv3 = domain.velocity_res.v3 as usize;
296
297        let dv = domain.dv();
298        let lv = [
299            domain.velocity.v1.to_f64().unwrap(),
300            domain.velocity.v2.to_f64().unwrap(),
301            domain.velocity.v3.to_f64().unwrap(),
302        ];
303
304        // Velocity spread: cold but resolved
305        let sigma = 0.3 * dv[0].max(dv[1]).max(dv[2]);
306        let norm = self.mean_density / ((2.0 * PI).sqrt() * sigma).powi(3);
307
308        // Strides for row-major 6D layout [x1, x2, x3, v1, v2, v3]
309        let s_v3 = 1usize;
310        let s_v2 = nv3;
311        let s_v1 = nv2 * nv3;
312        let s_x3 = nv1 * s_v1;
313        let s_x2 = nx3 * s_x3;
314        let s_x1 = nx2 * s_x2;
315
316        let total = nx1 * nx2 * nx3 * nv1 * nv2 * nv3;
317        let mut data = vec![0.0f64; total];
318
319        let inv_2sig2 = 1.0 / (2.0 * sigma * sigma);
320
321        for ix1 in 0..nx1 {
322            for ix2 in 0..nx2 {
323                for ix3 in 0..nx3 {
324                    let spatial_idx = ix1 * nx2 * nx3 + ix2 * nx3 + ix3;
325                    let v0x = vx0[spatial_idx];
326                    let v0y = vy0[spatial_idx];
327                    let v0z = vz0[spatial_idx];
328
329                    let base = ix1 * s_x1 + ix2 * s_x2 + ix3 * s_x3;
330
331                    for iv1 in 0..nv1 {
332                        let v1 = -lv[0] + (iv1 as f64 + 0.5) * dv[0];
333                        let dv1 = v1 - v0x;
334                        for iv2 in 0..nv2 {
335                            let v2 = -lv[1] + (iv2 as f64 + 0.5) * dv[1];
336                            let dv2 = v2 - v0y;
337                            for iv3 in 0..nv3 {
338                                let v3 = -lv[2] + (iv3 as f64 + 0.5) * dv[2];
339                                let dv3 = v3 - v0z;
340                                let v2sq = dv1 * dv1 + dv2 * dv2 + dv3 * dv3;
341                                let f = norm * (-v2sq * inv_2sig2).exp();
342                                data[base + iv1 * s_v1 + iv2 * s_v2 + iv3 * s_v3] = f;
343                            }
344                        }
345                    }
346                }
347            }
348        }
349
350        PhaseSpaceSnapshot {
351            data,
352            shape: [nx1, nx2, nx3, nv1, nv2, nv3],
353            time: 0.0,
354        }
355    }
356}
357
358/// 3D inverse FFT of a complex array, returning real part.
359/// Same pattern as `fft_3d_inverse` in `poisson/tensor_poisson.rs`.
360fn ifft_3d(data: &[rustfft::num_complex::Complex64], shape: [usize; 3]) -> Vec<f64> {
361    use rustfft::FftPlanner;
362    use rustfft::num_complex::Complex64;
363
364    let [n0, n1, n2] = shape;
365    let n_total = n0 * n1 * n2;
366    assert_eq!(data.len(), n_total);
367
368    let mut buf = data.to_vec();
369    let scale = 1.0 / n_total as f64;
370    let mut planner = FftPlanner::new();
371
372    // IFFT along axis 2 (fastest varying)
373    let ifft2 = planner.plan_fft_inverse(n2);
374    for i0 in 0..n0 {
375        for i1 in 0..n1 {
376            let start = i0 * n1 * n2 + i1 * n2;
377            ifft2.process(&mut buf[start..start + n2]);
378        }
379    }
380
381    // IFFT along axis 1
382    let ifft1 = planner.plan_fft_inverse(n1);
383    let mut line = vec![Complex64::new(0.0, 0.0); n1];
384    for i0 in 0..n0 {
385        for i2 in 0..n2 {
386            for i1 in 0..n1 {
387                line[i1] = buf[i0 * n1 * n2 + i1 * n2 + i2];
388            }
389            ifft1.process(&mut line);
390            for i1 in 0..n1 {
391                buf[i0 * n1 * n2 + i1 * n2 + i2] = line[i1];
392            }
393        }
394    }
395
396    // IFFT along axis 0
397    let ifft0 = planner.plan_fft_inverse(n0);
398    let mut line0 = vec![Complex64::new(0.0, 0.0); n0];
399    for i1 in 0..n1 {
400        for i2 in 0..n2 {
401            for i0 in 0..n0 {
402                line0[i0] = buf[i0 * n1 * n2 + i1 * n2 + i2];
403            }
404            ifft0.process(&mut line0);
405            for i0 in 0..n0 {
406                buf[i0 * n1 * n2 + i1 * n2 + i2] = line0[i0];
407            }
408        }
409    }
410
411    buf.iter().map(|c| c.re * scale).collect()
412}
413
414/// Simple single-mode Zel'dovich pancake IC for validation tests.
415pub struct ZeldovichSingleMode {
416    pub mean_density: f64,
417    pub amplitude: f64,
418    pub wavenumber: f64,
419    pub sigma_v: f64,
420}
421
422impl ZeldovichSingleMode {
423    /// Sample f(x,v) = ρ̄/(√(2π)σ_v)³ · exp(-(v - v₀(x))²/(2σ_v²))
424    /// where v₀(x) = (-A·sin(k·x₁), 0, 0)
425    pub fn sample_on_grid(&self, domain: &Domain) -> PhaseSpaceSnapshot {
426        use std::f64::consts::PI;
427
428        let nx1 = domain.spatial_res.x1 as usize;
429        let nx2 = domain.spatial_res.x2 as usize;
430        let nx3 = domain.spatial_res.x3 as usize;
431        let nv1 = domain.velocity_res.v1 as usize;
432        let nv2 = domain.velocity_res.v2 as usize;
433        let nv3 = domain.velocity_res.v3 as usize;
434
435        let dx = domain.dx();
436        let dv = domain.dv();
437        let lx = [
438            domain.spatial.x1.to_f64().unwrap(),
439            domain.spatial.x2.to_f64().unwrap(),
440            domain.spatial.x3.to_f64().unwrap(),
441        ];
442        let lv = [
443            domain.velocity.v1.to_f64().unwrap(),
444            domain.velocity.v2.to_f64().unwrap(),
445            domain.velocity.v3.to_f64().unwrap(),
446        ];
447
448        let s_v3 = 1usize;
449        let s_v2 = nv3;
450        let s_v1 = nv2 * nv3;
451        let s_x3 = nv1 * s_v1;
452        let s_x2 = nx3 * s_x3;
453        let s_x1 = nx2 * s_x2;
454
455        let total = nx1 * nx2 * nx3 * nv1 * nv2 * nv3;
456        let mut data = vec![0.0f64; total];
457
458        let sigma = self.sigma_v;
459        let norm = self.mean_density / ((2.0 * PI).sqrt() * sigma).powi(3);
460
461        for ix1 in 0..nx1 {
462            let x1 = -lx[0] + (ix1 as f64 + 0.5) * dx[0];
463            let v0x = -self.amplitude * (self.wavenumber * x1).sin();
464
465            for ix2 in 0..nx2 {
466                for ix3 in 0..nx3 {
467                    let base = ix1 * s_x1 + ix2 * s_x2 + ix3 * s_x3;
468
469                    for iv1 in 0..nv1 {
470                        let v1 = -lv[0] + (iv1 as f64 + 0.5) * dv[0];
471                        for iv2 in 0..nv2 {
472                            let v2 = -lv[1] + (iv2 as f64 + 0.5) * dv[1];
473                            for iv3 in 0..nv3 {
474                                let v3 = -lv[2] + (iv3 as f64 + 0.5) * dv[2];
475                                let dv1 = v1 - v0x;
476                                let dv2 = v2;
477                                let dv3 = v3;
478                                let v2sq = dv1 * dv1 + dv2 * dv2 + dv3 * dv3;
479                                let f = norm * (-v2sq / (2.0 * sigma * sigma)).exp();
480                                data[base + iv1 * s_v1 + iv2 * s_v2 + iv3 * s_v3] = f;
481                            }
482                        }
483                    }
484                }
485            }
486        }
487
488        PhaseSpaceSnapshot {
489            data,
490            shape: [nx1, nx2, nx3, nv1, nv2, nv3],
491            time: 0.0,
492        }
493    }
494}
495
496/// Growth factor D(a) for flat ΛCDM cosmology.
497///
498/// Uses the Carroll-Press-Turner (1992) approximation:
499///   D(a) ≈ (5/2) Ω_m(a) / [Ω_m(a)^{4/7} - Ω_Λ(a) + (1 + Ω_m(a)/2)(1 + Ω_Λ(a)/70)]
500///
501/// where Ω_m(a) = Ω_m / (Ω_m + Ω_Λ a³) and Ω_Λ(a) = Ω_Λ a³ / (Ω_m + Ω_Λ a³).
502/// Normalized so D(a) → a in the matter-dominated era.
503pub fn growth_factor(a: f64, params: &CosmologyParams) -> f64 {
504    let om = params.omega_m;
505    let ol = params.omega_lambda;
506
507    // E²(a) = Ω_m/a³ + Ω_Λ  (flat universe, Ω_k = 0)
508    let a3 = a * a * a;
509    let e2 = om / a3 + ol;
510
511    // Ω_m(a) and Ω_Λ(a) at scale factor a
512    let om_a = om / (a3 * e2);
513    let ol_a = ol / e2;
514
515    // Carroll-Press-Turner approximation
516    let d_unnorm =
517        2.5 * om_a / (om_a.powf(4.0 / 7.0) - ol_a + (1.0 + om_a / 2.0) * (1.0 + ol_a / 70.0));
518
519    // Normalize: in matter domination D ∝ a, so normalize to D(1) and scale
520    let e2_1 = om + ol;
521    let om_1 = om / e2_1;
522    let ol_1 = ol / e2_1;
523    let d_1 = 2.5 * om_1 / (om_1.powf(4.0 / 7.0) - ol_1 + (1.0 + om_1 / 2.0) * (1.0 + ol_1 / 70.0));
524
525    a * d_unnorm / d_1
526}
527
528#[cfg(test)]
529mod tests {
530    use super::*;
531    use crate::tooling::core::init::domain::{Domain, SpatialBoundType, VelocityBoundType};
532
533    #[test]
534    fn zeldovich_multi_mode_displacement_field() {
535        let domain = Domain::builder()
536            .spatial_extent(10.0)
537            .velocity_extent(2.0)
538            .spatial_resolution(8)
539            .velocity_resolution(8)
540            .t_final(1.0)
541            .spatial_bc(SpatialBoundType::Periodic)
542            .velocity_bc(VelocityBoundType::Open)
543            .build()
544            .unwrap();
545        let cosmo = CosmologyParams {
546            h0: 70.0,
547            omega_m: 0.3,
548            omega_lambda: 0.7,
549            a_init: 0.01,
550        };
551        let ic = ZeldovichIC::new(1.0, cosmo, 42);
552        let [sx, sy, sz] = ic.displacement_field(&domain);
553        let n = 8 * 8 * 8;
554        assert_eq!(sx.len(), n);
555        assert_eq!(sy.len(), n);
556        assert_eq!(sz.len(), n);
557        // Displacement should be finite and have zero mean (no DC mode)
558        assert!(sx.iter().all(|v| v.is_finite()));
559        let mean: f64 = sx.iter().sum::<f64>() / n as f64;
560        assert!(
561            mean.abs() < 1.0,
562            "displacement mean should be near zero, got {mean}"
563        );
564    }
565
566    #[test]
567    fn zeldovich_multi_mode_sample() {
568        // Use small H0 and late scale factor so velocities stay within grid extent.
569        // v ~ a * H(a) * f * s  where s ~ 0.1 * dx ~ 0.25, need v << velocity_extent.
570        let domain = Domain::builder()
571            .spatial_extent(10.0)
572            .velocity_extent(2.0)
573            .spatial_resolution(8)
574            .velocity_resolution(8)
575            .t_final(1.0)
576            .spatial_bc(SpatialBoundType::Periodic)
577            .velocity_bc(VelocityBoundType::Open)
578            .build()
579            .unwrap();
580        let cosmo = CosmologyParams {
581            h0: 1.0,
582            omega_m: 0.3,
583            omega_lambda: 0.7,
584            a_init: 1.0,
585        };
586        let ic = ZeldovichIC::new(1.0, cosmo, 42);
587        let snap = ic.sample_on_grid(&domain);
588        assert_eq!(snap.shape, [8, 8, 8, 8, 8, 8]);
589        assert!(snap.data.iter().all(|v| v.is_finite() && *v >= 0.0));
590        let mass: f64 = snap.data.iter().sum::<f64>();
591        assert!(mass > 0.0, "total mass should be positive");
592    }
593}