Skip to main content

sparse_ir/
basis.rs

1//! Finite temperature basis for SparseIR
2//!
3//! This module provides the `FiniteTempBasis` type which represents the
4//! intermediate representation (IR) basis for a given temperature.
5
6use std::sync::Arc;
7
8use crate::kernel::{CentrosymmKernel, KernelProperties, LogisticKernel};
9use crate::poly::{PiecewiseLegendrePolyVector, default_sampling_points};
10use crate::polyfourier::PiecewiseLegendreFTVector;
11use crate::sve::{SVEResult, TworkType, compute_sve};
12use crate::traits::{Bosonic, Fermionic, StatisticsType};
13
14// Re-export Statistics enum for C-API
15#[derive(Debug, Clone, Copy, PartialEq, Eq)]
16pub enum Statistics {
17    Fermionic,
18    Bosonic,
19}
20
21/// Finite temperature basis for imaginary time/frequency Green's functions
22///
23/// For a continuation kernel `K` from real frequencies `ω ∈ [-ωmax, ωmax]` to
24/// imaginary time `τ ∈ [0, β]`, this type stores the truncated singular
25/// value expansion or IR basis:
26///
27/// ```text
28/// K(τ, ω) ≈ sum(u[l](τ) * s[l] * v[l](ω) for l in 1:L)
29/// ```
30///
31/// This basis is inferred from a reduced form by appropriate scaling of
32/// the variables.
33///
34/// # Type Parameters
35///
36/// * `K` - Kernel type implementing `KernelProperties + CentrosymmKernel`
37/// * `S` - Statistics type (`Fermionic` or `Bosonic`)
38#[derive(Clone)]
39pub struct FiniteTempBasis<K, S>
40where
41    K: KernelProperties + CentrosymmKernel + Clone + 'static,
42    S: StatisticsType,
43{
44    /// The kernel used to construct this basis
45    kernel: K,
46
47    /// The SVE result (in scaled variables)
48    sve_result: Arc<SVEResult>,
49
50    /// Accuracy of the basis (relative error)
51    accuracy: f64,
52
53    /// Inverse temperature β
54    beta: f64,
55
56    /// Left singular functions on imaginary time axis τ ∈ [0, β]
57    /// Arc for efficient sharing (large immutable data)
58    u: Arc<PiecewiseLegendrePolyVector>,
59
60    /// Right singular functions on real frequency axis ω ∈ [-ωmax, ωmax]
61    /// Arc for efficient sharing (large immutable data)
62    v: Arc<PiecewiseLegendrePolyVector>,
63
64    /// Singular values
65    s: Vec<f64>,
66
67    /// Left singular functions on Matsubara frequency axis (Fourier transform of u)
68    /// Arc for efficient sharing (large immutable data)
69    uhat: Arc<PiecewiseLegendreFTVector<S>>,
70
71    /// Full uhat (before truncation to basis size)
72    /// Arc for efficient sharing (large immutable data, used for Matsubara sampling)
73    uhat_full: Arc<PiecewiseLegendreFTVector<S>>,
74
75    _phantom: std::marker::PhantomData<S>,
76}
77
78impl<K, S> FiniteTempBasis<K, S>
79where
80    K: KernelProperties + CentrosymmKernel + Clone + 'static,
81    S: StatisticsType,
82{
83    // ========== Getters ==========
84
85    /// Get a reference to the kernel
86    pub fn kernel(&self) -> &K {
87        &self.kernel
88    }
89
90    /// Get the SVE result
91    pub fn sve_result(&self) -> &Arc<SVEResult> {
92        &self.sve_result
93    }
94
95    /// Get the accuracy of the basis
96    pub fn accuracy(&self) -> f64 {
97        self.accuracy
98    }
99
100    /// Get the inverse temperature β
101    pub fn beta(&self) -> f64 {
102        self.beta
103    }
104
105    /// Get the left singular functions (u) on imaginary time axis
106    pub fn u(&self) -> &Arc<PiecewiseLegendrePolyVector> {
107        &self.u
108    }
109
110    /// Get the right singular functions (v) on real frequency axis
111    pub fn v(&self) -> &Arc<PiecewiseLegendrePolyVector> {
112        &self.v
113    }
114
115    /// Get the singular values
116    pub fn s(&self) -> &[f64] {
117        &self.s
118    }
119
120    /// Get the left singular functions on Matsubara frequency axis (uhat)
121    pub fn uhat(&self) -> &Arc<PiecewiseLegendreFTVector<S>> {
122        &self.uhat
123    }
124
125    /// Get the full uhat (before truncation)
126    pub fn uhat_full(&self) -> &Arc<PiecewiseLegendreFTVector<S>> {
127        &self.uhat_full
128    }
129
130    // ========== Other methods ==========
131
132    /// Get the frequency cutoff ωmax
133    pub fn wmax(&self) -> f64 {
134        self.kernel.lambda() / self.beta
135    }
136
137    /// Get default Matsubara sampling points as i64 indices (for C-API)
138    pub fn default_matsubara_sampling_points_i64(&self, positive_only: bool) -> Vec<i64>
139    where
140        S: 'static,
141    {
142        let freqs = self.default_matsubara_sampling_points(positive_only);
143        freqs.into_iter().map(|f| f.n()).collect()
144    }
145
146    /// Get default Matsubara sampling points as i64 indices with mitigate parameter (for C-API)
147    ///
148    /// # Panics
149    /// Panics if the kernel is not centrosymmetric. This method relies on
150    /// centrosymmetry to generate sampling points.
151    pub fn default_matsubara_sampling_points_i64_with_mitigate(
152        &self,
153        positive_only: bool,
154        mitigate: bool,
155        n_points: usize,
156    ) -> Vec<i64>
157    where
158        S: 'static,
159    {
160        if !self.kernel().is_centrosymmetric() {
161            panic!(
162                "default_matsubara_sampling_points_i64_with_mitigate is not supported for non-centrosymmetric kernels. \
163                 The current implementation relies on centrosymmetry to generate sampling points."
164            );
165        }
166        let fence = mitigate;
167        let freqs = Self::default_matsubara_sampling_points_impl(
168            &self.uhat_full,
169            n_points,
170            fence,
171            positive_only,
172        );
173        freqs.into_iter().map(|f| f.n()).collect()
174    }
175
176    /// Create a new FiniteTempBasis
177    ///
178    /// # Arguments
179    ///
180    /// * `kernel` - Kernel implementing `KernelProperties + CentrosymmKernel`
181    /// * `beta` - Inverse temperature (β > 0)
182    /// * `epsilon` - Accuracy parameter (optional, defaults to NaN for auto)
183    /// * `max_size` - Maximum number of basis functions (optional)
184    ///
185    /// # Returns
186    ///
187    /// A new FiniteTempBasis
188    pub fn new(kernel: K, beta: f64, epsilon: Option<f64>, max_size: Option<usize>) -> Self {
189        // Validate inputs
190        if beta <= 0.0 {
191            panic!("Inverse temperature beta must be positive, got {}", beta);
192        }
193
194        // Compute SVE
195        let epsilon_value = epsilon.unwrap_or(f64::NAN);
196        let sve_result = compute_sve(
197            kernel.clone(),
198            epsilon_value,
199            None, // cutoff
200            max_size,
201            TworkType::Auto,
202        );
203
204        Self::from_sve_result(kernel, beta, sve_result, epsilon, max_size)
205    }
206
207    /// Create basis from existing SVE result
208    ///
209    /// This is useful when you want to reuse the same SVE computation
210    /// for both fermionic and bosonic bases.
211    pub fn from_sve_result(
212        kernel: K,
213        beta: f64,
214        sve_result: SVEResult,
215        epsilon: Option<f64>,
216        max_size: Option<usize>,
217    ) -> Self {
218        // Get truncated u, s, v from SVE result
219        let (u_sve, s_sve, v_sve) = sve_result.part(epsilon, max_size);
220
221        // Calculate accuracy
222        let accuracy = if sve_result.s.len() > s_sve.len() {
223            sve_result.s[s_sve.len()] / sve_result.s[0]
224        } else {
225            sve_result.s[sve_result.s.len() - 1] / sve_result.s[0]
226        };
227
228        // Get kernel parameters
229        let lambda = kernel.lambda();
230        let omega_max = lambda / beta;
231
232        // Scale polynomials to new variables
233        // tau = β/2 * (x + 1), w = ωmax * y
234
235        // Transform u: x ∈ [-1, 1] → τ ∈ [0, β]
236        let u_knots: Vec<f64> = u_sve.get_polys()[0]
237            .knots
238            .iter()
239            .map(|&x| beta / 2.0 * (x + 1.0))
240            .collect();
241        let u_delta_x: Vec<f64> = u_sve.get_polys()[0]
242            .delta_x
243            .iter()
244            .map(|&dx| beta / 2.0 * dx)
245            .collect();
246        let u_symm: Vec<i32> = u_sve.get_polys().iter().map(|p| p.symm).collect();
247
248        let u = u_sve.rescale_domain(u_knots, Some(u_delta_x), Some(u_symm));
249
250        // Transform v: y ∈ [-1, 1] → ω ∈ [-ωmax, ωmax]
251        let v_knots: Vec<f64> = v_sve.get_polys()[0]
252            .knots
253            .iter()
254            .map(|&y| omega_max * y)
255            .collect();
256        let v_delta_x: Vec<f64> = v_sve.get_polys()[0]
257            .delta_x
258            .iter()
259            .map(|&dy| omega_max * dy)
260            .collect();
261        let v_symm: Vec<i32> = v_sve.get_polys().iter().map(|p| p.symm).collect();
262
263        let v = v_sve.rescale_domain(v_knots, Some(v_delta_x), Some(v_symm));
264
265        // Scale singular values
266        // s_scaled = sqrt(β/2 * ωmax) * ωmax^(-ypower) * s_sve
267        let ypower = kernel.ypower();
268        let scale_factor = (beta / 2.0 * omega_max).sqrt() * omega_max.powi(-ypower);
269        let s: Vec<f64> = s_sve.iter().map(|&x| scale_factor * x).collect();
270
271        // Construct uhat (Fourier transform of u)
272        // HACK: Fourier transforms only work on unit interval, so we scale the data
273        let uhat_base_full = sve_result.u.scale_data(beta.sqrt());
274        let conv_rad = kernel.conv_radius();
275
276        // Create statistics marker instance using Default trait
277        // S is a zero-sized type (ZST) like Fermionic or Bosonic
278        let stat_marker = S::default();
279
280        let uhat_full = PiecewiseLegendreFTVector::<S>::from_poly_vector(
281            &uhat_base_full,
282            stat_marker,
283            Some(conv_rad),
284        );
285
286        // Truncate uhat to basis size
287        let uhat_polyvec: Vec<_> = uhat_full.polyvec.iter().take(s.len()).cloned().collect();
288        let uhat = PiecewiseLegendreFTVector::from_vector(uhat_polyvec);
289
290        Self {
291            kernel,
292            sve_result: Arc::new(sve_result),
293            accuracy,
294            beta,
295            u: Arc::new(u),
296            v: Arc::new(v),
297            s,
298            uhat: Arc::new(uhat),
299            uhat_full: Arc::new(uhat_full),
300            _phantom: std::marker::PhantomData,
301        }
302    }
303
304    /// Get the size of the basis (number of basis functions)
305    pub fn size(&self) -> usize {
306        self.s.len()
307    }
308
309    /// Get the cutoff parameter Λ = β * ωmax
310    pub fn lambda(&self) -> f64 {
311        self.kernel.lambda()
312    }
313
314    /// Get the frequency cutoff ωmax
315    pub fn omega_max(&self) -> f64 {
316        self.lambda() / self.beta
317    }
318
319    /// Get significance of each singular value (s[i] / s[0])
320    pub fn significance(&self) -> Vec<f64> {
321        let s0 = self.s[0];
322        self.s.iter().map(|&s| s / s0).collect()
323    }
324
325    /// Get default tau sampling points
326    ///
327    /// Returns sampling points in imaginary time τ ∈ [-β/2, β/2].
328    ///
329    /// Roots are found with symmetry exploitation (matching Python 1.x / Julia v1),
330    /// then mapped to [-β/2, β/2] by folding τ_physical ∈ [0, β] around β/2.
331    pub fn default_tau_sampling_points(&self) -> Vec<f64> {
332        let points = self.default_tau_sampling_points_size_requested(self.size());
333        let basis_size = self.size();
334        if points.len() < basis_size {
335            eprintln!(
336                "Warning: Number of tau sampling points ({}) is less than basis size ({}). \
337                 Basis parameters: beta={}, wmax={}, epsilon={:.2e}",
338                points.len(),
339                basis_size,
340                self.beta,
341                self.wmax(),
342                self.accuracy()
343            );
344        }
345        points
346    }
347
348    /// Get default tau sampling points with a requested size
349    ///
350    /// Returns sampling points in τ ∈ [-β/2, β/2].
351    pub fn default_tau_sampling_points_size_requested(&self, size_requested: usize) -> Vec<f64> {
352        let x = default_sampling_points(&self.sve_result.u, size_requested);
353        let half_beta = self.beta / 2.0;
354        // Map roots to physical tau ∈ [0, β], then fold to [-β/2, β/2]
355        let mut smpl_taus: Vec<f64> = x
356            .iter()
357            .map(|&xi| {
358                let tau = half_beta * (xi + 1.0);
359                if tau <= half_beta {
360                    tau
361                } else {
362                    tau - self.beta
363                }
364            })
365            .collect();
366        smpl_taus.sort_by(|a, b| a.partial_cmp(b).unwrap());
367        smpl_taus
368    }
369
370    /// Get default Matsubara frequency sampling points
371    ///
372    /// Returns sampling points as MatsubaraFreq objects based on extrema
373    /// of the Matsubara basis functions (same algorithm as C++/Julia).
374    ///
375    /// # Arguments
376    /// * `positive_only` - If true, returns only non-negative frequencies
377    ///
378    /// # Returns
379    /// Vector of Matsubara frequency sampling points
380    ///
381    /// # Panics
382    /// Panics if the kernel is not centrosymmetric. This method relies on
383    /// centrosymmetry to generate sampling points.
384    pub fn default_matsubara_sampling_points(
385        &self,
386        positive_only: bool,
387    ) -> Vec<crate::freq::MatsubaraFreq<S>>
388    where
389        S: 'static,
390    {
391        if !self.kernel().is_centrosymmetric() {
392            panic!(
393                "default_matsubara_sampling_points is not supported for non-centrosymmetric kernels. \
394                 The current implementation relies on centrosymmetry to generate sampling points."
395            );
396        }
397        let fence = false;
398        let points = Self::default_matsubara_sampling_points_impl(
399            &self.uhat_full,
400            self.size(),
401            fence,
402            positive_only,
403        );
404        let basis_size = self.size();
405        // For positive_only=true, we need 2*n_sampling_points >= basis_size
406        // For positive_only=false, we need n_sampling_points >= basis_size
407        let effective_points = if positive_only {
408            2 * points.len()
409        } else {
410            points.len()
411        };
412        if effective_points < basis_size {
413            eprintln!(
414                "Warning: Number of Matsubara sampling points ({}{}) is less than basis size ({}). \
415                 Basis parameters: beta={}, wmax={}, epsilon={:.2e}",
416                points.len(),
417                if positive_only { " × 2" } else { "" },
418                basis_size,
419                self.beta,
420                self.wmax(),
421                self.accuracy()
422            );
423        }
424        points
425    }
426
427    /// Fence Matsubara sampling points to improve conditioning
428    ///
429    /// This function adds additional sampling points near the outer frequencies
430    /// to improve the conditioning of the sampling matrix. This is particularly
431    /// important for Matsubara sampling where we cannot freely choose sampling points.
432    ///
433    /// Implementation matches C++ version in `basis.hpp` (lines 407-452).
434    fn fence_matsubara_sampling(
435        omega_n: &mut Vec<crate::freq::MatsubaraFreq<S>>,
436        positive_only: bool,
437    ) where
438        S: StatisticsType + 'static,
439    {
440        use crate::freq::{BosonicFreq, MatsubaraFreq};
441
442        if omega_n.is_empty() {
443            return;
444        }
445
446        // Collect outer frequencies
447        let mut outer_frequencies = Vec::new();
448        if positive_only {
449            outer_frequencies.push(omega_n[omega_n.len() - 1]);
450        } else {
451            outer_frequencies.push(omega_n[0]);
452            outer_frequencies.push(omega_n[omega_n.len() - 1]);
453        }
454
455        for wn_outer in outer_frequencies {
456            let outer_val = wn_outer.n();
457            // In SparseIR.jl-v1, ωn_diff is always created as BosonicFreq
458            // This ensures diff_val is always even (valid for Bosonic)
459            let mut diff_val = 2 * (0.025 * outer_val as f64).round() as i64;
460
461            // Handle edge case: if diff_val is 0, set it to 2 (minimum even value for Bosonic)
462            if diff_val == 0 {
463                diff_val = 2;
464            }
465
466            // Get the n value from BosonicFreq (same as diff_val since it's even)
467            let wn_diff = BosonicFreq::new(diff_val).unwrap().n();
468
469            // Sign function: returns +1 if n > 0, -1 if n < 0, 0 if n == 0
470            // Matches C++ implementation: (a.get_n() > 0) - (a.get_n() < 0)
471            let sign_val = if outer_val > 0 {
472                1
473            } else if outer_val < 0 {
474                -1
475            } else {
476                0
477            };
478
479            // Check original size before adding (C++ checks wn.size() before each push)
480            let original_size = omega_n.len();
481            if original_size >= 20 {
482                // For Fermionic: wn_outer.n is odd, wn_diff is even, so wn_outer.n ± wn_diff is odd (valid)
483                // For Bosonic: wn_outer.n is even, wn_diff is even, so wn_outer.n ± wn_diff is even (valid)
484                let new_n = outer_val - sign_val * wn_diff;
485                if let Ok(new_freq) = MatsubaraFreq::<S>::new(new_n) {
486                    omega_n.push(new_freq);
487                }
488            }
489            if original_size >= 42 {
490                let new_n = outer_val + sign_val * wn_diff;
491                if let Ok(new_freq) = MatsubaraFreq::<S>::new(new_n) {
492                    omega_n.push(new_freq);
493                }
494            }
495        }
496
497        // Sort and remove duplicates using BTreeSet
498        let omega_n_set: std::collections::BTreeSet<MatsubaraFreq<S>> = omega_n.drain(..).collect();
499        *omega_n = omega_n_set.into_iter().collect();
500    }
501
502    pub fn default_matsubara_sampling_points_impl(
503        uhat_full: &PiecewiseLegendreFTVector<S>,
504        l: usize,
505        fence: bool,
506        positive_only: bool,
507    ) -> Vec<crate::freq::MatsubaraFreq<S>>
508    where
509        S: StatisticsType + 'static,
510    {
511        use crate::freq::MatsubaraFreq;
512        use crate::polyfourier::{find_extrema, sign_changes};
513        use std::collections::BTreeSet;
514
515        let mut l_requested = l;
516
517        // Adjust l_requested based on statistics (same as C++)
518        if S::STATISTICS == crate::traits::Statistics::Fermionic && l_requested % 2 != 0 {
519            l_requested += 1;
520        } else if S::STATISTICS == crate::traits::Statistics::Bosonic && l_requested % 2 == 0 {
521            l_requested += 1;
522        }
523
524        // Choose sign_changes or find_extrema based on l_requested
525        let mut omega_n = if l_requested < uhat_full.len() {
526            sign_changes(&uhat_full[l_requested], positive_only)
527        } else {
528            find_extrema(&uhat_full[uhat_full.len() - 1], positive_only)
529        };
530
531        // For bosons, include zero frequency explicitly to prevent conditioning issues
532        if S::STATISTICS == crate::traits::Statistics::Bosonic {
533            omega_n.push(MatsubaraFreq::<S>::new(0).unwrap());
534        }
535
536        // Sort and remove duplicates using BTreeSet
537        let omega_n_set: BTreeSet<MatsubaraFreq<S>> = omega_n.into_iter().collect();
538        let mut omega_n: Vec<MatsubaraFreq<S>> = omega_n_set.into_iter().collect();
539
540        // Check expected size
541        let expected_size = if positive_only {
542            l_requested.div_ceil(2)
543        } else {
544            l_requested
545        };
546
547        if omega_n.len() != expected_size {
548            eprintln!(
549                "Warning: Requested {} sampling frequencies for basis size L = {}, but got {}.",
550                expected_size,
551                l,
552                omega_n.len()
553            );
554        }
555
556        // Apply fencing if requested (same as C++ implementation)
557        if fence {
558            Self::fence_matsubara_sampling(&mut omega_n, positive_only);
559        }
560
561        omega_n
562    }
563    /// Get default omega (real frequency) sampling points
564    ///
565    /// Returns sampling points on the real-frequency axis ω ∈ [-ωmax, ωmax].
566    /// These are used as pole locations for the Discrete Lehmann Representation (DLR).
567    ///
568    /// The sampling points are chosen as the roots of the L-th basis function
569    /// in the spectral domain (v), which provides near-optimal conditioning.
570    ///
571    /// # Returns
572    /// Vector of real-frequency sampling points in [-ωmax, ωmax]
573    pub fn default_omega_sampling_points(&self) -> Vec<f64> {
574        let sz = self.size();
575
576        // Use UNTRUNCATED sve_result.v (same as C++)
577        // C++: default_sampling_points(*(sve_result->v), sz)
578        let y = default_sampling_points(&self.sve_result.v, sz);
579
580        // Scale to [-ωmax, ωmax]
581        let wmax = self.kernel.lambda() / self.beta;
582        let omega_points: Vec<f64> = y.into_iter().map(|yi| wmax * yi).collect();
583
584        omega_points
585    }
586}
587
588// ============================================================================
589// Trait implementations
590// ============================================================================
591
592impl<K, S> crate::basis_trait::Basis<S> for FiniteTempBasis<K, S>
593where
594    K: KernelProperties + CentrosymmKernel + Clone + 'static,
595    S: StatisticsType + 'static,
596{
597    type Kernel = K;
598
599    fn kernel(&self) -> &Self::Kernel {
600        &self.kernel
601    }
602
603    fn beta(&self) -> f64 {
604        self.beta
605    }
606
607    fn wmax(&self) -> f64 {
608        self.kernel.lambda() / self.beta
609    }
610
611    fn lambda(&self) -> f64 {
612        self.kernel.lambda()
613    }
614
615    fn size(&self) -> usize {
616        self.size()
617    }
618
619    fn accuracy(&self) -> f64 {
620        self.accuracy
621    }
622
623    fn significance(&self) -> Vec<f64> {
624        if let Some(&first_s) = self.s.first() {
625            self.s.iter().map(|&s| s / first_s).collect()
626        } else {
627            vec![]
628        }
629    }
630
631    fn svals(&self) -> Vec<f64> {
632        self.s.clone()
633    }
634
635    fn default_tau_sampling_points(&self) -> Vec<f64> {
636        self.default_tau_sampling_points()
637    }
638
639    fn default_matsubara_sampling_points(
640        &self,
641        positive_only: bool,
642    ) -> Vec<crate::freq::MatsubaraFreq<S>> {
643        self.default_matsubara_sampling_points(positive_only)
644    }
645
646    fn evaluate_tau(&self, tau: &[f64]) -> mdarray::DTensor<f64, 2> {
647        use crate::taufuncs::normalize_tau;
648        use mdarray::DTensor;
649
650        let n_points = tau.len();
651        let basis_size = self.size();
652
653        // Evaluate each basis function at all tau points
654        // Result: matrix[i, l] = u_l(tau[i])
655        // Note: tau can be in [-beta, beta] and will be normalized to [0, beta]
656        // self.u polynomials are already scaled to tau ∈ [0, beta] domain
657        DTensor::<f64, 2>::from_fn([n_points, basis_size], |idx| {
658            let i = idx[0]; // tau index
659            let l = idx[1]; // basis function index
660
661            // Normalize tau to [0, beta] with statistics-dependent sign
662            let (tau_norm, sign) = normalize_tau::<S>(tau[i], self.beta);
663
664            // Evaluate basis function directly (u polynomials are in tau domain)
665            sign * self.u[l].evaluate(tau_norm)
666        })
667    }
668
669    fn evaluate_matsubara(
670        &self,
671        freqs: &[crate::freq::MatsubaraFreq<S>],
672    ) -> mdarray::DTensor<num_complex::Complex<f64>, 2> {
673        use mdarray::DTensor;
674        use num_complex::Complex;
675
676        let n_points = freqs.len();
677        let basis_size = self.size();
678
679        // Evaluate each basis function at all Matsubara frequencies
680        // Result: matrix[i, l] = uhat_l(iωn[i])
681        DTensor::<Complex<f64>, 2>::from_fn([n_points, basis_size], |idx| {
682            let i = idx[0]; // frequency index
683            let l = idx[1]; // basis function index
684            self.uhat[l].evaluate(&freqs[i])
685        })
686    }
687
688    fn evaluate_omega(&self, omega: &[f64]) -> mdarray::DTensor<f64, 2> {
689        use mdarray::DTensor;
690
691        let n_points = omega.len();
692        let basis_size = self.size();
693
694        // Evaluate each spectral basis function at all omega points
695        // Result: matrix[i, l] = V_l(omega[i])
696        DTensor::<f64, 2>::from_fn([n_points, basis_size], |idx| {
697            let i = idx[0]; // omega index
698            let l = idx[1]; // basis function index
699            self.v[l].evaluate(omega[i])
700        })
701    }
702
703    fn default_omega_sampling_points(&self) -> Vec<f64> {
704        self.default_omega_sampling_points()
705    }
706}
707
708// ============================================================================
709// Type aliases
710// ============================================================================
711
712/// Type alias for fermionic basis with LogisticKernel
713pub type FermionicBasis = FiniteTempBasis<LogisticKernel, Fermionic>;
714
715/// Type alias for bosonic basis with LogisticKernel
716pub type BosonicBasis = FiniteTempBasis<LogisticKernel, Bosonic>;
717
718#[cfg(test)]
719#[path = "basis_tests.rs"]
720mod basis_tests;