1use num_complex::Complex64;
7use std::f64::consts::PI;
8
9use crate::freq::MatsubaraFreq;
10use crate::poly::{PiecewiseLegendrePoly, PiecewiseLegendrePolyVector};
11use crate::special_functions::spherical_bessel_j;
12use crate::traits::{Bosonic, Fermionic, Statistics, StatisticsType};
13
14#[derive(Debug, Clone)]
16pub struct PowerModel {
17 pub moments: Vec<f64>,
18}
19
20impl PowerModel {
21 pub fn new(moments: Vec<f64>) -> Self {
23 Self { moments }
24 }
25}
26
27#[derive(Debug, Clone)]
32pub struct PiecewiseLegendreFT<S: StatisticsType> {
33 pub poly: PiecewiseLegendrePoly,
35 pub n_asymp: f64,
37 pub model: PowerModel,
39 _phantom: std::marker::PhantomData<S>,
40}
41
42pub type FermionicPiecewiseLegendreFT = PiecewiseLegendreFT<Fermionic>;
44pub type BosonicPiecewiseLegendreFT = PiecewiseLegendreFT<Bosonic>;
45
46impl<S: StatisticsType> PiecewiseLegendreFT<S> {
47 pub fn new(poly: PiecewiseLegendrePoly, _stat: S, n_asymp: Option<f64>) -> Self {
57 if (poly.xmin - (-1.0)).abs() > 1e-12 || (poly.xmax - 1.0).abs() > 1e-12 {
59 panic!("Only interval [-1, 1] is supported for Fourier transform");
60 }
61
62 let n_asymp = n_asymp.unwrap_or(f64::INFINITY);
63 let model = Self::power_model(&poly);
64
65 Self {
66 poly,
67 n_asymp,
68 model,
69 _phantom: std::marker::PhantomData,
70 }
71 }
72
73 pub fn get_n_asymp(&self) -> f64 {
75 self.n_asymp
76 }
77
78 pub fn get_statistics(&self) -> Statistics {
80 S::STATISTICS
81 }
82
83 pub fn zeta(&self) -> i64 {
85 match S::STATISTICS {
86 Statistics::Fermionic => 1,
87 Statistics::Bosonic => 0,
88 }
89 }
90
91 pub fn get_poly(&self) -> &PiecewiseLegendrePoly {
93 &self.poly
94 }
95
96 pub fn evaluate(&self, omega: &MatsubaraFreq<S>) -> Complex64 {
104 let n = omega.get_n() as i32;
105 if (n as f64).abs() < self.n_asymp {
106 self.compute_unl_inner(&self.poly, n)
107 } else {
108 self.giw(n)
109 }
110 }
111
112 pub fn evaluate_at_n(&self, n: i64) -> Complex64 {
114 match MatsubaraFreq::<S>::new(n) {
115 Ok(omega) => self.evaluate(&omega),
116 Err(_) => Complex64::new(0.0, 0.0), }
118 }
119
120 pub fn evaluate_at_ns(&self, ns: &[i64]) -> Vec<Complex64> {
122 ns.iter().map(|&n| self.evaluate_at_n(n)).collect()
123 }
124
125 fn power_model(poly: &PiecewiseLegendrePoly) -> PowerModel {
127 let deriv_x1 = poly.derivs(1.0);
128 let moments = Self::power_moments(&deriv_x1, poly.l);
129 PowerModel::new(moments)
130 }
131
132 fn power_moments(deriv_x1: &[f64], l: i32) -> Vec<f64> {
134 let statsign = match S::STATISTICS {
135 Statistics::Fermionic => -1.0,
136 Statistics::Bosonic => 1.0,
137 };
138
139 let mut moments = deriv_x1.to_vec();
140 for (m, moment) in moments.iter_mut().enumerate() {
141 let m_f64 = (m + 1) as f64; *moment *=
143 -(statsign * (-1.0_f64).powi(m_f64 as i32) + (-1.0_f64).powi(l)) / 2.0_f64.sqrt();
144 }
145 moments
146 }
147
148 fn compute_unl_inner(&self, poly: &PiecewiseLegendrePoly, wn: i32) -> Complex64 {
150 let wred = PI / 4.0 * wn as f64;
151 let phase_wi = Self::phase_stable(poly, wn);
152 let mut res = Complex64::new(0.0, 0.0);
153
154 let order_max = poly.polyorder;
155 let segment_count = poly.knots.len() - 1;
156
157 for order in 0..order_max {
158 for j in 0..segment_count {
159 let data_oj = poly.data[[order, j]];
160 let tnl = Self::get_tnl(order as i32, wred * poly.delta_x[j]);
161 res += data_oj * tnl * phase_wi[j] / poly.norms[j];
162 }
163 }
164
165 res / 2.0_f64.sqrt()
166 }
167
168 fn giw(&self, wn: i32) -> Complex64 {
170 let iw = Complex64::new(0.0, PI / 2.0 * wn as f64);
171 if wn == 0 {
172 return Complex64::new(0.0, 0.0);
173 }
174
175 let inv_iw = 1.0 / iw;
176
177 inv_iw * Self::evalpoly(inv_iw, &self.model.moments)
178 }
179
180 fn evalpoly(x: Complex64, coeffs: &[f64]) -> Complex64 {
182 let mut result = Complex64::new(0.0, 0.0);
183 for i in (0..coeffs.len()).rev() {
184 result = result * x + Complex64::new(coeffs[i], 0.0);
185 }
186 result
187 }
188
189 fn shift_xmid(knots: &[f64], delta_x: &[f64]) -> (Vec<f64>, Vec<i32>) {
195 let n_segments = delta_x.len();
196 let delta_x_half: Vec<f64> = delta_x.iter().map(|&dx| dx / 2.0).collect();
197
198 let mut xmid_m1 = Vec::with_capacity(n_segments);
200 let mut cumsum = 0.0;
201 for i in 0..n_segments {
202 cumsum += delta_x[i];
203 xmid_m1.push(cumsum - delta_x_half[i]);
204 }
205
206 let mut xmid_p1 = Vec::with_capacity(n_segments);
208 let mut cumsum_rev = 0.0;
209 for i in (0..n_segments).rev() {
210 cumsum_rev += delta_x[i];
211 xmid_p1.insert(0, -cumsum_rev + delta_x_half[i]);
212 }
213
214 let xmid_0: Vec<f64> = (0..n_segments)
216 .map(|i| knots[i + 1] - delta_x_half[i])
217 .collect();
218
219 let mut xmid_diff = Vec::with_capacity(n_segments);
221 let mut extra_shift = Vec::with_capacity(n_segments);
222
223 for i in 0..n_segments {
224 let shift = xmid_0[i].round() as i32;
225 extra_shift.push(shift);
226
227 let diff = match shift {
229 -1 => xmid_m1[i],
230 0 => xmid_0[i],
231 1 => xmid_p1[i],
232 _ => xmid_0[i], };
234 xmid_diff.push(diff);
235 }
236
237 (xmid_diff, extra_shift)
238 }
239
240 fn phase_stable(poly: &PiecewiseLegendrePoly, wn: i32) -> Vec<Complex64> {
245 let (xmid_diff, extra_shift) = Self::shift_xmid(&poly.knots, &poly.delta_x);
246 let mut phase_wi = Vec::with_capacity(xmid_diff.len());
247
248 let im_unit = Complex64::new(0.0, 1.0);
249
250 for j in 0..xmid_diff.len() {
251 let power = ((wn * (extra_shift[j] + 1)) % 4 + 4) % 4; let im_power = im_unit.powi(power);
254
255 let arg = PI * wn as f64 * xmid_diff[j] / 2.0;
257 let cispi = Complex64::new(arg.cos(), arg.sin());
258
259 phase_wi.push(im_power * cispi);
260 }
261
262 phase_wi
263 }
264
265 pub fn sign_changes(&self, positive_only: bool) -> Vec<MatsubaraFreq<S>> {
273 let f = Self::func_for_part(self);
274 let x0 = Self::find_all_roots(&f, DEFAULT_GRID);
275
276 let mut matsubara_indices: Vec<i64> = x0.into_iter().map(|x| 2 * x + self.zeta()).collect();
278
279 if !positive_only {
280 Self::symmetrize_matsubara_inplace(&mut matsubara_indices);
281 }
282
283 matsubara_indices
284 .into_iter()
285 .filter_map(|n| MatsubaraFreq::<S>::new(n).ok())
286 .collect()
287 }
288
289 pub fn find_extrema(&self, positive_only: bool) -> Vec<MatsubaraFreq<S>> {
297 let f = Self::func_for_part(self);
298 let x0 = Self::discrete_extrema(&f, DEFAULT_GRID);
299
300 let mut matsubara_indices: Vec<i64> = x0.into_iter().map(|x| 2 * x + self.zeta()).collect();
302
303 if !positive_only {
304 Self::symmetrize_matsubara_inplace(&mut matsubara_indices);
305 }
306
307 matsubara_indices
308 .into_iter()
309 .filter_map(|n| MatsubaraFreq::<S>::new(n).ok())
310 .collect()
311 }
312
313 fn func_for_part(&self) -> impl Fn(i64) -> f64 + '_ {
315 let parity = self.poly.symm;
316 let poly_ft = self.clone();
317
318 move |n| {
319 let omega = match MatsubaraFreq::<S>::new(n) {
320 Ok(omega) => omega,
321 Err(_) => return 0.0,
322 };
323 let value = poly_ft.evaluate(&omega);
324
325 let result = match parity {
326 1 => match S::STATISTICS {
327 Statistics::Fermionic => value.im,
328 Statistics::Bosonic => value.re,
329 },
330 -1 => match S::STATISTICS {
331 Statistics::Fermionic => value.re,
332 Statistics::Bosonic => value.im,
333 },
334 0 => {
335 value.re
337 }
338 _ => panic!("Cannot detect parity for symm = {}", parity),
339 };
340
341 if n == 0 || n == 1 || n == 2 {
343 println!("n={}, value={}, result={}", n, value, result);
344 }
345
346 result
347 }
348 }
349
350 fn find_all_roots<F>(f: &F, xgrid: &[i64]) -> Vec<i64>
352 where
353 F: Fn(i64) -> f64,
354 {
355 if xgrid.is_empty() {
356 return Vec::new();
357 }
358
359 let fx: Vec<f64> = xgrid.iter().map(|&x| f(x)).collect();
361
362 let mut x_hit = Vec::new();
364 for i in 0..fx.len() {
365 if fx[i] == 0.0 {
366 x_hit.push(xgrid[i]);
367 }
368 }
369
370 let mut sign_change = Vec::new();
372 for i in 0..fx.len() - 1 {
373 let has_sign_change = fx[i].signum() != fx[i + 1].signum();
374 let not_hit = fx[i] != 0.0 && fx[i + 1] != 0.0;
375 let both_nonzero = fx[i].abs() > 1e-12 && fx[i + 1].abs() > 1e-12;
376 sign_change.push(has_sign_change && not_hit && both_nonzero);
377 }
378
379 if sign_change.iter().all(|&sc| !sc) {
381 x_hit.sort();
382 return x_hit;
383 }
384
385 let mut a_intervals = Vec::new();
387 let mut b_intervals = Vec::new();
388 let mut fa_values = Vec::new();
389
390 for i in 0..sign_change.len() {
391 if sign_change[i] {
392 a_intervals.push(xgrid[i]);
393 b_intervals.push(xgrid[i + 1]);
394 fa_values.push(fx[i]);
395 }
396 }
397
398 for i in 0..a_intervals.len() {
400 let root = Self::bisect(&f, a_intervals[i], b_intervals[i], fa_values[i]);
401 x_hit.push(root);
402 }
403
404 x_hit.sort();
406 x_hit
407 }
408
409 fn bisect<F>(f: &F, a: i64, b: i64, fa: f64) -> i64
411 where
412 F: Fn(i64) -> f64,
413 {
414 let mut a = a;
415 let mut b = b;
416 let mut fa = fa;
417
418 loop {
419 if (b - a).abs() <= 1 {
420 return a;
421 }
422
423 let mid = (a + b) / 2;
424 let fmid = f(mid);
425
426 if fa.signum() != fmid.signum() {
427 b = mid;
428 } else {
429 a = mid;
430 fa = fmid;
431 }
432 }
433 }
434
435 fn discrete_extrema<F>(f: &F, xgrid: &[i64]) -> Vec<i64>
437 where
438 F: Fn(i64) -> f64,
439 {
440 if xgrid.len() < 3 {
441 return Vec::new();
442 }
443
444 let fx: Vec<f64> = xgrid.iter().map(|&x| f(x)).collect();
445 let mut extrema = Vec::new();
446
447 for i in 1..fx.len() - 1 {
449 let prev = fx[i - 1];
450 let curr = fx[i];
451 let next = fx[i + 1];
452
453 if curr > prev && curr > next {
455 extrema.push(xgrid[i]);
456 }
457 else if curr < prev && curr < next {
459 extrema.push(xgrid[i]);
460 }
461 }
462
463 extrema
464 }
465
466 fn symmetrize_matsubara_inplace(xs: &mut Vec<i64>) {
468 xs.retain(|&x| x != 0);
470
471 xs.sort();
473
474 let negatives: Vec<i64> = xs.iter().rev().map(|&x| -x).collect();
476
477 xs.splice(0..0, negatives);
479 }
480
481 pub fn get_tnl(l: i32, w: f64) -> Complex64 {
487 let abs_w = w.abs();
488
489 let sph_bessel = spherical_bessel_j(l, abs_w);
491
492 let im_unit = Complex64::new(0.0, 1.0);
494 let im_power = im_unit.powi(l);
495 let result = 2.0 * im_power * sph_bessel;
496
497 if w < 0.0 { result.conj() } else { result }
499 }
500}
501
502#[derive(Debug, Clone)]
504pub struct PiecewiseLegendreFTVector<S: StatisticsType> {
505 pub polyvec: Vec<PiecewiseLegendreFT<S>>,
506 _phantom: std::marker::PhantomData<S>,
507}
508
509pub type FermionicPiecewiseLegendreFTVector = PiecewiseLegendreFTVector<Fermionic>;
511pub type BosonicPiecewiseLegendreFTVector = PiecewiseLegendreFTVector<Bosonic>;
512
513impl<S: StatisticsType> PiecewiseLegendreFTVector<S> {
514 pub fn new() -> Self {
516 Self {
517 polyvec: Vec::new(),
518 _phantom: std::marker::PhantomData,
519 }
520 }
521
522 pub fn from_vector(polyvec: Vec<PiecewiseLegendreFT<S>>) -> Self {
524 Self {
525 polyvec,
526 _phantom: std::marker::PhantomData,
527 }
528 }
529
530 pub fn len(&self) -> usize {
532 self.polyvec.len()
533 }
534
535 pub fn is_empty(&self) -> bool {
537 self.polyvec.is_empty()
538 }
539
540 pub fn from_poly_vector(
542 polys: &PiecewiseLegendrePolyVector,
543 _stat: S,
544 n_asymp: Option<f64>,
545 ) -> Self {
546 let mut polyvec = Vec::with_capacity(polys.size());
547
548 for i in 0..polys.size() {
549 let poly = polys.get(i).unwrap().clone();
550 let ft_poly = PiecewiseLegendreFT::new(poly, _stat, n_asymp);
551 polyvec.push(ft_poly);
552 }
553
554 Self {
555 polyvec,
556 _phantom: std::marker::PhantomData,
557 }
558 }
559
560 pub fn size(&self) -> usize {
562 self.polyvec.len()
563 }
564
565 pub fn get(&self, index: usize) -> Option<&PiecewiseLegendreFT<S>> {
567 self.polyvec.get(index)
568 }
569
570 pub fn get_mut(&mut self, index: usize) -> Option<&mut PiecewiseLegendreFT<S>> {
572 self.polyvec.get_mut(index)
573 }
574
575 pub fn set(&mut self, index: usize, poly: PiecewiseLegendreFT<S>) -> Result<(), String> {
577 if index >= self.polyvec.len() {
578 return Err(format!("Index {} out of range", index));
579 }
580 self.polyvec[index] = poly;
581 Ok(())
582 }
583
584 pub fn similar(&self) -> Self {
586 Self::new()
587 }
588
589 pub fn n_asymp(&self) -> f64 {
591 self.polyvec.first().map_or(f64::INFINITY, |p| p.n_asymp)
592 }
593
594 pub fn evaluate_at(&self, omega: &MatsubaraFreq<S>) -> Vec<Complex64> {
596 self.polyvec
597 .iter()
598 .map(|poly| poly.evaluate(omega))
599 .collect()
600 }
601
602 pub fn evaluate_at_many(&self, omegas: &[MatsubaraFreq<S>]) -> Vec<Vec<Complex64>> {
604 omegas.iter().map(|omega| self.evaluate_at(omega)).collect()
605 }
606}
607
608impl<S: StatisticsType> std::ops::Index<usize> for PiecewiseLegendreFTVector<S> {
610 type Output = PiecewiseLegendreFT<S>;
611
612 fn index(&self, index: usize) -> &Self::Output {
613 &self.polyvec[index]
614 }
615}
616
617impl<S: StatisticsType> std::ops::IndexMut<usize> for PiecewiseLegendreFTVector<S> {
618 fn index_mut(&mut self, index: usize) -> &mut Self::Output {
619 &mut self.polyvec[index]
620 }
621}
622
623impl<S: StatisticsType> Default for PiecewiseLegendreFTVector<S> {
625 fn default() -> Self {
626 Self::new()
627 }
628}
629
630const DEFAULT_GRID: &[i64] = &[
636 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 25,
637 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49,
638 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63, 64, 65, 66, 68, 69, 71, 72, 74, 76, 77,
639 79, 81, 82, 84, 86, 88, 90, 92, 94, 96, 98, 100, 103, 105, 107, 109, 112, 114, 117, 119, 122,
640 125, 128, 130, 133, 136, 139, 142, 145, 148, 152, 155, 158, 162, 165, 169, 173, 177, 181, 184,
641 189, 193, 197, 201, 206, 210, 215, 219, 224, 229, 234, 239, 245, 250, 256, 261, 267, 273, 279,
642 285, 291, 297, 304, 311, 317, 324, 331, 339, 346, 354, 362, 369, 378, 386, 394, 403, 412, 421,
643 430, 439, 449, 459, 469, 479, 490, 501, 512, 523, 534, 546, 558, 570, 583, 595, 608, 622, 635,
644 649, 663, 678, 693, 708, 724, 739, 756, 772, 789, 806, 824, 842, 861, 879, 899, 918, 939, 959,
645 980, 1002, 1024, 1046, 1069, 1092, 1116, 1141, 1166, 1191, 1217, 1244, 1271, 1299, 1327, 1357,
646 1386, 1417, 1448, 1479, 1512, 1545, 1579, 1613, 1649, 1685, 1722, 1759, 1798, 1837, 1878, 1919,
647 1961, 2004, 2048, 2092, 2138, 2185, 2233, 2282, 2332, 2383, 2435, 2488, 2543, 2599, 2655, 2714,
648 2773, 2834, 2896, 2959, 3024, 3090, 3158, 3227, 3298, 3370, 3444, 3519, 3596, 3675, 3756, 3838,
649 3922, 4008, 4096, 4185, 4277, 4371, 4466, 4564, 4664, 4766, 4870, 4977, 5086, 5198, 5311, 5428,
650 5547, 5668, 5792, 5919, 6049, 6181, 6316, 6455, 6596, 6741, 6888, 7039, 7193, 7351, 7512, 7676,
651 7844, 8016, 8192, 8371, 8554, 8742, 8933, 9129, 9328, 9533, 9741, 9955, 10173, 10396, 10623,
652 10856, 11094, 11336, 11585, 11838, 12098, 12363, 12633, 12910, 13193, 13482, 13777, 14078,
653 14387, 14702, 15024, 15353, 15689, 16032, 16384, 16742, 17109, 17484, 17866, 18258, 18657,
654 19066, 19483, 19910, 20346, 20792, 21247, 21712, 22188, 22673, 23170, 23677, 24196, 24726,
655 25267, 25820, 26386, 26964, 27554, 28157, 28774, 29404, 30048, 30706, 31378, 32065, 32768,
656 33485, 34218, 34968, 35733, 36516, 37315, 38132, 38967, 39821, 40693, 41584, 42494, 43425,
657 44376, 45347, 46340, 47355, 48392, 49452, 50535, 51641, 52772, 53928, 55108, 56315, 57548,
658 58809, 60096, 61412, 62757, 64131, 65536, 66971, 68437, 69936, 71467, 73032, 74631, 76265,
659 77935, 79642, 81386, 83168, 84989, 86850, 88752, 90695, 92681, 94711, 96785, 98904, 101070,
660 103283, 105545, 107856, 110217, 112631, 115097, 117618, 120193, 122825, 125514, 128263, 131072,
661 133942, 136875, 139872, 142935, 146064, 149263, 152531, 155871, 159284, 162772, 166337, 169979,
662 173701, 177504, 181391, 185363, 189422, 193570, 197809, 202140, 206566, 211090, 215712, 220435,
663 225262, 230195, 235236, 240387, 245650, 251029, 256526, 262144, 267884, 273750, 279744, 285870,
664 292129, 298526, 305063, 311743, 318569, 325545, 332674, 339958, 347402, 355009, 362783, 370727,
665 378845, 387141, 395618, 404281, 413133, 422180, 431424, 440871, 450525, 460390, 470472, 480774,
666 491301, 502059, 513053, 524288, 535768, 547500, 559488, 571740, 584259, 597053, 610126, 623487,
667 637139, 651091, 665348, 679917, 694805, 710019, 725567, 741455, 757690, 774282, 791236, 808562,
668 826267, 844360, 862849, 881743, 901051, 920781, 940944, 961548, 982603, 1004119, 1026107,
669 1048576, 1071536, 1095000, 1118977, 1143480, 1168519, 1194106, 1220253, 1246974, 1274279,
670 1302182, 1330696, 1359834, 1389611, 1420039, 1451134, 1482910, 1515381, 1548564, 1582473,
671 1617125, 1652535, 1688721, 1725699, 1763487, 1802102, 1841563, 1881888, 1923096, 1965207,
672 2008239, 2052214, 2097152, 2143073, 2190000, 2237955, 2286960, 2337038, 2388212, 2440507,
673 2493948, 2548558, 2604364, 2661392, 2719669, 2779222, 2840079, 2902269, 2965820, 3030763,
674 3097128, 3164947, 3234250, 3305071, 3377443, 3451399, 3526975, 3604205, 3683127, 3763777,
675 3846193, 3930414, 4016479, 4104428, 4194304, 4286147, 4380001, 4475911, 4573920, 4674076,
676 4776425, 4881015, 4987896, 5097116, 5208729, 5322785, 5439339, 5558445, 5680159, 5804538,
677 5931641, 6061527, 6194257, 6329894, 6468501, 6610142, 6754886, 6902798, 7053950, 7208411,
678 7366255, 7527555, 7692387, 7860828, 8032958, 8208857, 8388608, 8572294, 8760003, 8951822,
679 9147841, 9348153, 9552851, 9762031, 9975792, 10194233, 10417458, 10645571, 10878678, 11116890,
680 11360318, 11609077, 11863283, 12123055, 12388515, 12659788, 12937002, 13220285, 13509772,
681 13805597, 14107900, 14416823, 14732510, 15055110, 15384774, 15721657, 16065917, 16417714,
682 16777216, 17144589, 17520006, 17903645, 18295683, 18696307, 19105702, 19524063, 19951584,
683 20388467, 20834916, 21291142, 21757357, 22233781, 22720637, 23218155, 23726566, 24246110,
684 24777031, 25319577, 25874004, 26440571, 27019544, 27611195, 28215801, 28833647, 29465021,
685 30110221, 30769549, 31443315, 32131834, 32835429, 33554432,
686];
687
688pub fn sign_changes<S: StatisticsType + 'static>(
692 u_hat: &PiecewiseLegendreFT<S>,
693 positive_only: bool,
694) -> Vec<MatsubaraFreq<S>> {
695 let f = func_for_part(u_hat);
696 let x0 = find_all(&f, DEFAULT_GRID);
697
698 let mut indices: Vec<i64> = x0.iter().map(|&x| 2 * x + u_hat.zeta()).collect();
700
701 if !positive_only {
702 symmetrize_matsubara_inplace(&mut indices);
703 }
704
705 indices
706 .iter()
707 .filter_map(|&n| MatsubaraFreq::<S>::new(n).ok())
708 .collect()
709}
710
711pub fn find_extrema<S: StatisticsType + 'static>(
715 u_hat: &PiecewiseLegendreFT<S>,
716 positive_only: bool,
717) -> Vec<MatsubaraFreq<S>> {
718 let f = func_for_part(u_hat);
719 let x0 = discrete_extrema(&f, DEFAULT_GRID);
720
721 let mut indices: Vec<i64> = x0.iter().map(|&x| 2 * x + u_hat.zeta()).collect();
723
724 if !positive_only {
725 symmetrize_matsubara_inplace(&mut indices);
726 }
727
728 indices
729 .iter()
730 .filter_map(|&n| MatsubaraFreq::<S>::new(n).ok())
731 .collect()
732}
733
734fn func_for_part<S: StatisticsType + 'static>(
736 poly_ft: &PiecewiseLegendreFT<S>,
737) -> Box<dyn Fn(i64) -> f64> {
738 let parity = poly_ft.poly.symm();
739 let zeta = poly_ft.zeta();
740
741 let poly_ft_clone = poly_ft.clone();
743
744 Box::new(move |n: i64| {
745 let omega = MatsubaraFreq::<S>::new(2 * n + zeta).unwrap();
746 let value = poly_ft_clone.evaluate(&omega);
747
748 if parity == 1 {
750 if S::STATISTICS == Statistics::Bosonic {
752 value.re
753 } else {
754 value.im
755 }
756 } else if parity == -1 {
757 if S::STATISTICS == Statistics::Bosonic {
759 value.im
760 } else {
761 value.re
762 }
763 } else {
764 panic!("Cannot detect parity");
765 }
766 })
767}
768
769fn find_all(f: &dyn Fn(i64) -> f64, xgrid: &[i64]) -> Vec<i64> {
771 if xgrid.is_empty() {
772 return Vec::new();
773 }
774
775 let mut results = Vec::new();
776 let mut prev_val = f(xgrid[0]);
777
778 for i in 1..xgrid.len() {
779 let val = f(xgrid[i]);
780 if prev_val.signum() != val.signum() && val != 0.0 {
781 results.push(xgrid[i - 1]);
782 }
783 prev_val = val;
784 }
785
786 results
787}
788
789fn discrete_extrema(f: &dyn Fn(i64) -> f64, xgrid: &[i64]) -> Vec<i64> {
791 let mut results = Vec::new();
792
793 if xgrid.len() < 3 {
794 return results;
795 }
796
797 let mut prev_val = f(xgrid[0]);
798 let mut curr_val = f(xgrid[1]);
799
800 for i in 2..xgrid.len() {
801 let next_val = f(xgrid[i]);
802
803 if (curr_val > prev_val && curr_val > next_val)
805 || (curr_val < prev_val && curr_val < next_val)
806 {
807 results.push(xgrid[i - 1]);
808 }
809
810 prev_val = curr_val;
811 curr_val = next_val;
812 }
813
814 results
815}
816
817fn symmetrize_matsubara_inplace(xs: &mut Vec<i64>) {
819 if xs.is_empty() {
820 return;
821 }
822
823 xs.retain(|&x| x != 0);
825
826 let positives: Vec<i64> = xs.iter().filter(|&&x| x > 0).copied().collect();
828 let mut negatives: Vec<i64> = positives.iter().map(|&x| -x).collect();
829
830 xs.append(&mut negatives);
831 xs.sort();
832 xs.dedup();
833}
834
835#[cfg(test)]
836#[path = "polyfourier_tests.rs"]
837mod polyfourier_tests;