sparse_ir/kernel.rs
1//! Kernel implementations for SparseIR
2//!
3//! This module provides kernel implementations for analytical continuation in quantum many-body physics.
4//! The kernels are used in Fredholm integral equations of the first kind.
5//!
6//! u(x) = integral of K(x, y) v(y) dy
7//!
8//! where x ∈ [xmin, xmax] and y ∈ [ymin, ymax].
9//!
10//! In general, the kernel is applied to a scaled spectral function rho'(y) as:
11//!
12//! integral of K(x, y) rho'(y) dy,
13//!
14//! where ρ'(y) = w(y) ρ(y). The weight function w(y) transforms the original spectral
15//! function ρ(y) into the scaled version ρ'(y) used in the integral equation.
16
17use crate::numeric::CustomNumeric;
18use crate::traits::{Statistics, StatisticsType};
19use std::fmt::Debug;
20
21#[derive(Debug, Clone, Copy, PartialEq, Eq)]
22pub enum SymmetryType {
23 Even,
24 Odd,
25}
26
27impl SymmetryType {
28 pub fn sign(self) -> i32 {
29 match self {
30 SymmetryType::Even => 1,
31 SymmetryType::Odd => -1,
32 }
33 }
34}
35
36/// Trait for SVE (Singular Value Expansion) hints
37///
38/// Provides discretization hints for singular value expansion of a given kernel.
39/// This includes segment information and numerical parameters for efficient computation.
40pub trait SVEHints<T>: Debug + Send + Sync
41where
42 T: Copy + Debug + Send + Sync,
43{
44 /// Get the x-axis segments for discretization
45 ///
46 /// For centrosymmetric kernels, returns only positive values (x >= 0) including the endpoints.
47 /// The returned vector contains segments from [0, xmax] where xmax is the
48 /// upper bound of the x domain.
49 ///
50 /// For non-centrosymmetric kernels, returns segments covering the full domain
51 /// [-xmax, xmax].
52 fn segments_x(&self) -> Vec<T>;
53
54 /// Get the y-axis segments for discretization
55 ///
56 /// For centrosymmetric kernels, returns only positive values (y >= 0) including the endpoints.
57 /// The returned vector contains segments from [0, ymax] where ymax is the
58 /// upper bound of the y domain.
59 ///
60 /// For non-centrosymmetric kernels, returns segments covering the full domain
61 /// [-ymax, ymax].
62 fn segments_y(&self) -> Vec<T>;
63
64 /// Get the number of singular values hint
65 fn nsvals(&self) -> usize;
66
67 /// Get the number of Gauss points for quadrature
68 fn ngauss(&self) -> usize;
69}
70
71/// Trait for kernel type properties (static characteristics)
72pub trait KernelProperties {
73 /// Associated type for SVE hints
74 type SVEHintsType<T>: SVEHints<T> + Clone
75 where
76 T: Copy + Debug + Send + Sync + CustomNumeric + 'static;
77 /// Power with which the y coordinate scales.
78 ///
79 /// For most kernels, this is 0 (no scaling).
80 /// For RegularizedBoseKernel, this is 1 (linear scaling).
81 fn ypower(&self) -> i32;
82
83 /// Convergence radius of the Matsubara basis asymptotic model
84 ///
85 /// For improved numerical accuracy, IR basis functions on Matsubara axis
86 /// can be evaluated from asymptotic expression for |n| > conv_radius.
87 fn conv_radius(&self) -> f64;
88
89 /// Get the upper bound of the x domain
90 fn xmax(&self) -> f64;
91
92 /// Get the upper bound of the y domain
93 fn ymax(&self) -> f64;
94
95 /// A regularizer for a bosonic kernel for avoiding a divergence at omega = 0.
96 ///
97 /// The bosonic kernel diverges at omega = 0.
98 /// This function returns a regularizer w(beta, omega) that avoids this divergence.
99 ///
100 /// G(τ) = - ∫ K(τ, ω) ρ(ω) dω = - ∫ [K(τ, ω) w(ω)] [ρ(ω)/w(ω)] dω.
101 ///
102 /// The spectral function ρ(ω) and the weight function w(ω) must vanish linearly at omega = 0.
103 /// For a fermionic kernel, this function is expected to return 1.0.
104 ///
105 /// # Arguments
106 ///
107 /// * `beta` - Inverse temperature
108 /// * `omega` - Frequency
109 ///
110 /// # Returns
111 ///
112 /// The regularizer value w(beta, omega)
113 fn regularizer<S: StatisticsType + 'static>(&self, beta: f64, omega: f64) -> f64;
114
115 /// Create SVE hints for this kernel type.
116 ///
117 /// Provides discretization hints for singular value expansion computation.
118 /// The hints include segment information and numerical parameters optimized
119 /// for the specific kernel type.
120 ///
121 /// @param epsilon Target accuracy for the SVE computation
122 /// @return SVE hints specific to this kernel type
123 fn sve_hints<T>(&self, epsilon: f64) -> Self::SVEHintsType<T>
124 where
125 T: Copy + Debug + Send + Sync + CustomNumeric + 'static;
126}
127
128/// Trait for general kernels (both centrosymmetric and non-centrosymmetric)
129///
130/// This trait provides the basic interface for computing kernel values.
131/// Centrosymmetric kernels should implement `CentrosymmKernel` instead,
132/// which provides additional optimizations.
133pub trait AbstractKernel: Send + Sync {
134 /// Compute the kernel value K(x, y) with high precision
135 ///
136 /// # Arguments
137 ///
138 /// * `x` - The x coordinate (typically in [-xmax, xmax])
139 /// * `y` - The y coordinate (typically in [-ymax, ymax])
140 ///
141 /// # Returns
142 ///
143 /// The kernel value K(x, y)
144 fn compute<T: CustomNumeric + Copy + Debug>(&self, x: T, y: T) -> T;
145
146 /// Check if the kernel is centrosymmetric
147 ///
148 /// Returns true if and only if K(x, y) == K(-x, -y) for all values of x and y.
149 /// This allows the kernel to be block-diagonalized, speeding up the
150 /// singular value expansion by a factor of 4.
151 ///
152 /// # Returns
153 ///
154 /// True if the kernel is centrosymmetric, false otherwise.
155 fn is_centrosymmetric(&self) -> bool {
156 false
157 }
158}
159
160/// Trait for centrosymmetric kernels
161///
162/// Centrosymmetric kernels satisfy K(x, y) = K(-x, -y) and can be decomposed
163/// into even and odd components for efficient computation.
164pub trait CentrosymmKernel: AbstractKernel {
165 /// Compute the reduced kernel value
166 ///
167 /// K_red(x, y) = K(x, y) + sign * K(x, -y)
168 /// where sign = 1 for even symmetry and sign = -1 for odd symmetry
169 fn compute_reduced<T: CustomNumeric + Copy + Debug>(
170 &self,
171 x: T,
172 y: T,
173 symmetry: SymmetryType,
174 ) -> T;
175
176 /// Get the cutoff parameter Λ (lambda)
177 fn lambda(&self) -> f64;
178}
179
180/// Logistic kernel for fermionic analytical continuation
181///
182/// This kernel implements K(x, y) = exp(-Λy(x + 1)/2)/(1 + exp(-Λy))
183/// where x ∈ [-1, 1] and y ∈ [-1, 1]
184#[derive(Debug, Clone, Copy)]
185pub struct LogisticKernel {
186 lambda: f64,
187}
188
189impl LogisticKernel {
190 /// Create a new logistic kernel with the given cutoff parameter
191 pub fn new(lambda: f64) -> Self {
192 Self { lambda }
193 }
194
195 /// Get the cutoff parameter
196 pub fn lambda(&self) -> f64 {
197 self.lambda
198 }
199}
200
201impl KernelProperties for LogisticKernel {
202 type SVEHintsType<T>
203 = LogisticSVEHints<T>
204 where
205 T: Copy + Debug + Send + Sync + CustomNumeric + 'static;
206 fn ypower(&self) -> i32 {
207 0 // No y-power scaling for LogisticKernel
208 }
209
210 fn conv_radius(&self) -> f64 {
211 40.0 * self.lambda
212 }
213
214 fn xmax(&self) -> f64 {
215 1.0
216 }
217 fn ymax(&self) -> f64 {
218 1.0
219 }
220
221 fn regularizer<S: StatisticsType + 'static>(&self, beta: f64, omega: f64) -> f64 {
222 match S::STATISTICS {
223 Statistics::Fermionic => {
224 // For fermionic statistics: regularizer = 1.0 (safe, no division by zero)
225 1.0
226 }
227 Statistics::Bosonic => {
228 // For bosonic statistics: regularizer = tanh(0.5 * beta * omega) (safe, handles omega=0 case)
229 // This avoids division by zero when tanh(0.5 * beta * omega) approaches zero
230 (0.5 * beta * omega).tanh()
231 }
232 }
233 }
234
235 fn sve_hints<T>(&self, epsilon: f64) -> Self::SVEHintsType<T>
236 where
237 T: Copy + Debug + Send + Sync + CustomNumeric + 'static,
238 {
239 LogisticSVEHints::new(*self, epsilon)
240 }
241}
242
243pub fn compute_logistic_kernel<T: CustomNumeric>(lambda: f64, x: T, y: T) -> T {
244 let x_plus: T = T::from_f64_unchecked(1.0) + x;
245 let x_minus: T = T::from_f64_unchecked(1.0) - x;
246
247 let u_plus: T = T::from_f64_unchecked(0.5) * x_plus;
248 let u_minus: T = T::from_f64_unchecked(0.5) * x_minus;
249 let v: T = T::from_f64_unchecked(lambda) * y;
250
251 let mabs_v: T = -v.abs_as_same_type();
252 let numerator: T = if v >= T::from_f64_unchecked(0.0) {
253 (u_plus * mabs_v).exp()
254 } else {
255 (u_minus * mabs_v).exp()
256 };
257 let denominator: T = T::from_f64_unchecked(1.0) + mabs_v.exp();
258 numerator / denominator
259}
260
261fn compute_logistic_kernel_reduced_odd<T: CustomNumeric>(lambda: f64, x: T, y: T) -> T {
262 // For x * y around 0, antisymmetrization introduces cancellation, which
263 // reduces the relative precision. To combat this, we replace the
264 // values with the explicit form
265 let v_half: T = T::from_f64_unchecked(lambda * 0.5) * y;
266 // Use direct comparison to match C++ implementation (avoid precision loss from to_f64())
267 let xy_small: bool = x * v_half < T::from_f64_unchecked(1.0);
268 let cosh_finite: bool = v_half < T::from_f64_unchecked(85.0);
269 if xy_small && cosh_finite {
270 -(v_half * x).sinh() / v_half.cosh()
271 } else {
272 let k_plus = compute_logistic_kernel(lambda, x, y);
273 let k_minus = compute_logistic_kernel(lambda, x, -y);
274 k_plus - k_minus
275 }
276}
277
278impl AbstractKernel for LogisticKernel {
279 fn compute<T: CustomNumeric + Copy + Debug>(&self, x: T, y: T) -> T {
280 compute_logistic_kernel(self.lambda, x, y)
281 }
282
283 fn is_centrosymmetric(&self) -> bool {
284 true
285 }
286}
287
288impl CentrosymmKernel for LogisticKernel {
289 fn compute_reduced<T: CustomNumeric + Copy + Debug>(
290 &self,
291 x: T,
292 y: T,
293 symmetry: SymmetryType,
294 ) -> T {
295 match symmetry {
296 SymmetryType::Even => self.compute(x, y) + self.compute(x, -y),
297 SymmetryType::Odd => compute_logistic_kernel_reduced_odd(self.lambda, x, y),
298 }
299 }
300
301 fn lambda(&self) -> f64 {
302 self.lambda
303 }
304}
305
306/// SVE hints for LogisticKernel
307#[derive(Debug, Clone)]
308pub struct LogisticSVEHints<T> {
309 kernel: LogisticKernel,
310 epsilon: f64,
311 _phantom: std::marker::PhantomData<T>,
312}
313
314impl<T> LogisticSVEHints<T>
315where
316 T: Copy + Debug + Send + Sync,
317{
318 pub fn new(kernel: LogisticKernel, epsilon: f64) -> Self {
319 Self {
320 kernel,
321 epsilon,
322 _phantom: std::marker::PhantomData,
323 }
324 }
325}
326
327impl<T> SVEHints<T> for LogisticSVEHints<T>
328where
329 T: Copy + Debug + Send + Sync + CustomNumeric,
330{
331 fn segments_x(&self) -> Vec<T> {
332 // Direct implementation that generates only non-negative sample points
333 // This is equivalent to the C++ implementation but without the full symmetric array creation
334 let nzeros = std::cmp::max((15.0 * self.kernel.lambda().log10()).round() as usize, 1);
335
336 // Create a range of values
337 let mut temp = vec![0.0; nzeros];
338 for i in 0..nzeros {
339 temp[i] = 0.143 * i as f64;
340 }
341
342 // Calculate diffs using the inverse hyperbolic cosine
343 let mut diffs = vec![0.0; nzeros];
344 for i in 0..nzeros {
345 diffs[i] = 1.0 / temp[i].cosh();
346 }
347
348 // Calculate cumulative sum of diffs
349 let mut zeros = vec![0.0; nzeros];
350 zeros[0] = diffs[0];
351 for i in 1..nzeros {
352 zeros[i] = zeros[i - 1] + diffs[i];
353 }
354
355 // Normalize zeros
356 let last_zero = zeros[nzeros - 1];
357 for i in 0..nzeros {
358 zeros[i] /= last_zero;
359 }
360
361 // Create segments with only non-negative values (x >= 0) including endpoints [0, xmax]
362 let mut segments = Vec::with_capacity(nzeros + 1);
363
364 // Add 0.0 endpoint
365 segments.push(<T as CustomNumeric>::from_f64_unchecked(0.0));
366
367 // Add positive zeros (already in [0, 1] range)
368 for i in 0..nzeros {
369 segments.push(<T as CustomNumeric>::from_f64_unchecked(zeros[i]));
370 }
371
372 // Ensure segments are sorted in ascending order [0, ..., xmax]
373 segments.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
374
375 segments
376 }
377
378 fn segments_y(&self) -> Vec<T> {
379 // Direct implementation that generates only non-negative sample points
380 // This is equivalent to the C++ implementation but without the full symmetric array creation
381 let nzeros = std::cmp::max((20.0 * self.kernel.lambda().log10()).round() as usize, 2);
382
383 // Initial differences (from C++ implementation)
384 let mut diffs = vec![
385 0.01523, 0.03314, 0.04848, 0.05987, 0.06703, 0.07028, 0.07030, 0.06791, 0.06391,
386 0.05896, 0.05358, 0.04814, 0.04288, 0.03795, 0.03342, 0.02932, 0.02565, 0.02239,
387 0.01951, 0.01699,
388 ];
389
390 // Truncate diffs if necessary
391 if nzeros < diffs.len() {
392 diffs.truncate(nzeros);
393 }
394
395 // Calculate trailing differences
396 for i in 20..nzeros {
397 let x = 0.141 * i as f64;
398 diffs.push(0.25 * (-x).exp());
399 }
400
401 // Calculate cumulative sum of diffs
402 let mut zeros = Vec::with_capacity(nzeros);
403 zeros.push(diffs[0]);
404 for i in 1..nzeros {
405 zeros.push(zeros[i - 1] + diffs[i]);
406 }
407
408 // Normalize zeros
409 let last_zero = zeros[nzeros - 1];
410 for i in 0..nzeros {
411 zeros[i] /= last_zero;
412 }
413 zeros.pop(); // Remove last element
414
415 // Updated nzeros
416 let nzeros = zeros.len();
417
418 // Adjust zeros
419 for i in 0..nzeros {
420 zeros[i] -= 1.0;
421 }
422
423 // Generate segments directly from negative zeros
424 let mut segments: Vec<T> = Vec::new();
425
426 segments.push(<T as CustomNumeric>::from_f64_unchecked(1.0));
427
428 // Add absolute values of negative zeros
429 for i in 0..nzeros {
430 let abs_val = -zeros[i];
431 segments.push(<T as CustomNumeric>::from_f64_unchecked(abs_val));
432 }
433
434 if segments[segments.len() - 1].abs_as_same_type() > T::epsilon() {
435 segments.push(<T as CustomNumeric>::from_f64_unchecked(0.0));
436 }
437
438 // Sort in ascending order
439 segments.sort_by(|a, b| a.partial_cmp(b).unwrap());
440
441 segments
442 }
443
444 fn nsvals(&self) -> usize {
445 let log10_lambda = self.kernel.lambda().log10().max(1.0);
446 ((25.0 + log10_lambda) * log10_lambda).round() as usize
447 }
448
449 fn ngauss(&self) -> usize {
450 if self.epsilon >= 1e-8 { 10 } else { 16 }
451 }
452}
453
454// ============================================================================
455// RegularizedBoseKernel
456// ============================================================================
457
458/// Regularized bosonic analytical continuation kernel
459///
460/// In dimensionless variables x = 2τ/β - 1, y = βω/Λ, the integral kernel is:
461///
462/// ```text
463/// K(x, y) = y * exp(-Λ y (x + 1) / 2) / (1 - exp(-Λ y))
464/// ```
465///
466/// This kernel is used for bosonic Green's functions. The factor y regularizes
467/// the singularity at ω = 0, making the kernel well-behaved for numerical work.
468///
469/// The dimensionalized kernel is related by:
470/// ```text
471/// K(τ, ω) = ωmax * K(2τ/β - 1, ω/ωmax)
472/// ```
473/// where ωmax = Λ/β.
474///
475/// # Properties
476/// - **Centrosymmetric**: K(x, y) = K(-x, -y)
477/// - **ypower = 1**: Spectral function transforms as ρ'(y) = y * ρ(y)
478/// - **Bosonic only**: Does not support fermionic statistics
479/// - **Weight function**: w(β, ω) = 1/ω for bosonic statistics
480///
481/// # Numerical Stability
482/// The expression v / (exp(v) - 1) is evaluated using expm1 for small |v|.
483#[derive(Debug, Clone, Copy, PartialEq)]
484pub struct RegularizedBoseKernel {
485 /// Kernel cutoff parameter Λ = β × ωmax
486 pub lambda: f64,
487}
488
489impl RegularizedBoseKernel {
490 /// Create a new RegularizedBoseKernel
491 ///
492 /// # Arguments
493 /// * `lambda` - Kernel cutoff Λ (must be non-negative)
494 ///
495 /// # Panics
496 /// Panics if lambda < 0
497 pub fn new(lambda: f64) -> Self {
498 if lambda < 0.0 {
499 panic!("Kernel cutoff Λ must be non-negative, got {}", lambda);
500 }
501 Self { lambda }
502 }
503
504 /// Compute kernel value with numerical stability
505 ///
506 /// Evaluates K(x, y) = y * exp(-Λy(x+1)/2) / (1 - exp(-Λy))
507 /// using expm1 for better accuracy near y = 0.
508 ///
509 /// # Arguments
510 /// * `x` - Normalized time coordinate (x ∈ [-1, 1])
511 /// * `y` - Normalized frequency coordinate (y ∈ [-1, 1])
512 /// * `x_plus` - Precomputed x + 1 (optional, for efficiency)
513 /// * `x_minus` - Precomputed 1 - x (optional, for efficiency)
514 fn compute_impl<T>(&self, x: T, y: T, x_plus: Option<T>, x_minus: Option<T>) -> T
515 where
516 T: CustomNumeric,
517 {
518 // Convert lambda and constants to type T at the beginning
519 let lambda_t = T::from_f64_unchecked(self.lambda);
520 let half = T::from_f64_unchecked(0.5);
521 let inv_lambda = T::from_f64_unchecked(1.0) / lambda_t;
522
523 // u_plus = (x + 1) / 2, u_minus = (1 - x) / 2
524 // x_plus and x_minus are (x+1) and (1-x), so we need to multiply by 0.5
525 let u_plus = x_plus
526 .map(|xp| half * xp)
527 .unwrap_or_else(|| half * (T::from_f64_unchecked(1.0) + x));
528 let u_minus = x_minus
529 .map(|xm| half * xm)
530 .unwrap_or_else(|| half * (T::from_f64_unchecked(1.0) - x));
531
532 let v = lambda_t * y;
533 let absv = v.abs_as_same_type();
534
535 // Handle y ≈ 0 using Taylor expansion
536 // K(x,y) = 1/Λ - xy/2 + (1/24)Λ(3x² - 1)y² + O(y³)
537 // The limit as y->0 is 1/Λ (using L'Hopital's rule)
538 // For |Λy| < 2e-14, use first-order approximation
539 // This avoids division by zero when exp(-|Λy|) ≈ 1
540 if absv.to_f64() < 2e-14 {
541 let term0 = inv_lambda;
542 let term1 = half * x * y;
543 return term0 - term1;
544 }
545
546 // enum_val = exp(-|v| * (v >= 0 ? u_plus : u_minus))
547 let enum_val = if v >= T::from_f64_unchecked(0.0) {
548 (-absv * u_plus).exp()
549 } else {
550 (-absv * u_minus).exp()
551 };
552
553 // Handle v / (exp(v) - 1) with numerical stability using expm1
554 // Follows SparseIR.jl implementation: denom = absv / expm1(-absv)
555 // This is more accurate than exp(-absv) - 1 for small arguments
556 let denom = if absv.to_f64() >= 1e-200 {
557 let expm1_neg_absv = CustomNumeric::exp_m1(-absv);
558 absv / expm1_neg_absv
559 } else {
560 // For very small absv, use -1 (matches SparseIR.jl: -one(absv))
561 -T::from_f64_unchecked(1.0)
562 };
563
564 // K(x, y) = -1/Λ * enum_val * denom
565 // Since denom is negative (exp(-absv) < 1), final result is positive
566 -inv_lambda * enum_val * denom
567 }
568}
569
570impl KernelProperties for RegularizedBoseKernel {
571 type SVEHintsType<T>
572 = RegularizedBoseSVEHints<T>
573 where
574 T: Copy + Debug + Send + Sync + CustomNumeric + 'static;
575
576 fn ypower(&self) -> i32 {
577 1 // Spectral function transforms as ρ'(y) = y * ρ(y)
578 }
579
580 fn conv_radius(&self) -> f64 {
581 40.0 * self.lambda
582 }
583
584 fn xmax(&self) -> f64 {
585 1.0
586 }
587
588 fn ymax(&self) -> f64 {
589 1.0
590 }
591
592 fn regularizer<S: StatisticsType + 'static>(&self, _beta: f64, omega: f64) -> f64 {
593 match S::STATISTICS {
594 Statistics::Fermionic => {
595 panic!("RegularizedBoseKernel does not support fermionic functions");
596 }
597 Statistics::Bosonic => {
598 // regularizer = ω (safe, no division)
599 omega
600 }
601 }
602 }
603
604 fn sve_hints<T>(&self, epsilon: f64) -> Self::SVEHintsType<T>
605 where
606 T: Copy + Debug + Send + Sync + CustomNumeric + 'static,
607 {
608 RegularizedBoseSVEHints::new(*self, epsilon)
609 }
610}
611
612impl AbstractKernel for RegularizedBoseKernel {
613 fn compute<T: CustomNumeric + Copy + Debug>(&self, x: T, y: T) -> T {
614 let x_plus = Some(T::from_f64_unchecked(1.0) + x);
615 let x_minus = Some(T::from_f64_unchecked(1.0) - x);
616 self.compute_impl(x, y, x_plus, x_minus)
617 }
618
619 fn is_centrosymmetric(&self) -> bool {
620 true
621 }
622}
623
624impl CentrosymmKernel for RegularizedBoseKernel {
625 fn compute_reduced<T: CustomNumeric + Copy + Debug>(
626 &self,
627 x: T,
628 y: T,
629 symmetry: SymmetryType,
630 ) -> T {
631 match symmetry {
632 SymmetryType::Even => self.compute(x, y) + self.compute(x, -y),
633 SymmetryType::Odd => {
634 // For RegularizedBoseKernel, use sinh formulation for numerical stability
635 // K(x,y) - K(x,-y) = -y * sinh(Λ y x / 2) / sinh(Λ y / 2)
636 let v_half = T::from_f64_unchecked(self.lambda * 0.5) * y;
637 let xv_half = x * v_half;
638 let xy_small = xv_half.to_f64().abs() < 1.0;
639 let sinh_finite = v_half.to_f64().abs() < 85.0 && v_half.to_f64().abs() > 1e-200;
640
641 if xy_small && sinh_finite {
642 // Use sinh formulation for numerical stability
643 // NOTE: Minus sign is critical! (matches Julia/C++ implementation)
644 -y * xv_half.sinh() / v_half.sinh()
645 } else {
646 // Fall back to direct computation
647 self.compute(x, y) - self.compute(x, -y)
648 }
649 }
650 }
651 }
652
653 fn lambda(&self) -> f64 {
654 self.lambda
655 }
656}
657
658/// SVE hints for RegularizedBoseKernel
659#[derive(Debug, Clone)]
660pub struct RegularizedBoseSVEHints<T> {
661 kernel: RegularizedBoseKernel,
662 epsilon: f64,
663 _phantom: std::marker::PhantomData<T>,
664}
665
666impl<T> RegularizedBoseSVEHints<T>
667where
668 T: Copy + Debug + Send + Sync,
669{
670 pub fn new(kernel: RegularizedBoseKernel, epsilon: f64) -> Self {
671 Self {
672 kernel,
673 epsilon,
674 _phantom: std::marker::PhantomData,
675 }
676 }
677}
678
679impl<T> SVEHints<T> for RegularizedBoseSVEHints<T>
680where
681 T: Copy + Debug + Send + Sync + CustomNumeric + 'static,
682{
683 fn segments_x(&self) -> Vec<T> {
684 // Returns segments for x >= 0 domain only
685 // C++/Julia: nzeros = max(round(15 * log10(lambda)), 15)
686 let nzeros = ((15.0 * self.kernel.lambda.log10()).round() as usize).max(15);
687
688 // temp[i] = 0.18 * i
689 let mut temp = vec![0.0; nzeros];
690 for i in 0..nzeros {
691 temp[i] = 0.18 * i as f64;
692 }
693
694 // diffs[i] = 1.0 / cosh(temp[i])
695 let mut diffs = vec![0.0; nzeros];
696 for i in 0..nzeros {
697 diffs[i] = 1.0 / temp[i].cosh();
698 }
699
700 // Cumulative sum
701 let mut zeros = vec![0.0; nzeros];
702 zeros[0] = diffs[0];
703 for i in 1..nzeros {
704 zeros[i] = zeros[i - 1] + diffs[i];
705 }
706
707 // Normalize
708 let last_zero = zeros[nzeros - 1];
709 for i in 0..nzeros {
710 zeros[i] /= last_zero;
711 }
712
713 // Create segments with only non-negative values: [0, zeros[0], zeros[1], ...]
714 let mut segments = Vec::with_capacity(nzeros + 1);
715 segments.push(T::from_f64_unchecked(0.0));
716 for i in 0..nzeros {
717 segments.push(T::from_f64_unchecked(zeros[i]));
718 }
719
720 // Ensure sorted (should already be sorted, but verify)
721 segments.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
722
723 segments
724 }
725
726 fn segments_y(&self) -> Vec<T> {
727 // Returns segments for y >= 0 domain only
728 // C++/Julia: nzeros = max(round(20 * log10(lambda)), 20)
729 let nzeros = ((20.0 * self.kernel.lambda.log10()).round() as usize).max(20);
730
731 // diffs[j] = 0.12 / exp(0.0337 * j * log(j + 1))
732 let mut diffs = vec![0.0; nzeros];
733 for j in 0..nzeros {
734 let j_f64 = j as f64;
735 let exponent = 0.0337 * j_f64 * (j_f64 + 1.0).ln();
736 diffs[j] = 0.12 / exponent.exp();
737 }
738
739 // Cumulative sum
740 let mut zeros = vec![0.0; nzeros];
741 zeros[0] = diffs[0];
742 for i in 1..nzeros {
743 zeros[i] = zeros[i - 1] + diffs[i];
744 }
745
746 // Normalize by last value, then remove it
747 let last_zero = zeros[nzeros - 1];
748 for i in 0..nzeros {
749 zeros[i] /= last_zero;
750 }
751 zeros.pop();
752
753 // After normalization and pop, zeros are in [0, 1) range
754 // Create segments: [0, zeros[0], zeros[1], ..., 1]
755 let mut segments = Vec::with_capacity(zeros.len() + 2);
756 segments.push(T::from_f64_unchecked(0.0));
757 for z in zeros {
758 segments.push(T::from_f64_unchecked(z));
759 }
760 segments.push(T::from_f64_unchecked(1.0));
761
762 // Ensure sorted (should already be sorted, but verify)
763 segments.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
764
765 segments
766 }
767
768 fn nsvals(&self) -> usize {
769 // C++: int(round(28 * max(1.0, log10(lambda))))
770 let log10_lambda = self.kernel.lambda.log10().max(1.0);
771 (28.0 * log10_lambda).round() as usize
772 }
773
774 fn ngauss(&self) -> usize {
775 if self.epsilon >= 1e-8 { 10 } else { 16 }
776 }
777}
778
779#[cfg(test)]
780#[path = "kernel_tests.rs"]
781mod tests;