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