ginger/rootfinding.rs
1use super::{Matrix2, Vector2};
2
3type Vec2 = Vector2<f64>;
4type Mat2 = Matrix2<f64>;
5
6const PI: f64 = std::f64::consts::PI;
7
8/// The below code defines a struct named Options with three fields: max_iters, tolerance, and tol_ind.
9///
10/// Properties:
11///
12/// * `max_iters`: The `max_iters` property represents the maximum number of iterations allowed for a
13/// certain algorithm or process. It is of type `usize`, which means it can only hold non-negative
14/// integer values.
15/// * `tolerance`: The `tolerance` property is a floating-point number that represents the tolerance for convergence
16/// in an algorithm. It is used to determine when the algorithm has reached a satisfactory solution.
17/// * `tol_ind`: The `tol_ind` property in the `Options` struct represents the tolerance for individual
18/// values. It is a floating-point number (`f64`) that determines the acceptable difference between the
19/// expected value and the actual value for each element in a calculation or comparison.
20#[derive(Debug)]
21pub struct Options {
22 pub max_iters: usize,
23 pub tolerance: f64,
24 pub tol_ind: f64,
25}
26
27/// The below code is implementing the `Default` trait for the `Options` struct in Rust. The `Default`
28/// trait provides a default value for a type, which can be used when creating an instance of the type
29/// without specifying any values. In this case, the `default` function is defined to return an instance
30/// of the `Options` struct with default values for the `max_iters`, `tolerance`, and `tol_ind` fields.
31impl Default for Options {
32 fn default() -> Self {
33 Options {
34 max_iters: 2000,
35 tolerance: 1e-12,
36 tol_ind: 1e-15,
37 }
38 }
39}
40
41/// The function `make_adjoint` calculates the adjoint matrix between two vectors.
42///
43/// Arguments:
44///
45/// * `vr`: A vector representing the direction of the reference frame's x-axis.
46/// * `vp`: The parameter `vp` represents a vector `vp = (p, s)`, where `p` and `s` are the components
47/// of the vector.
48///
49/// Returns:
50///
51/// The function `make_adjoint` returns a `Mat2` object.
52#[inline]
53pub fn make_adjoint(vr: &Vec2, vp: &Vec2) -> Mat2 {
54 let (r, q) = (vr.x_, vr.y_);
55 let (p, s) = (vp.x_, vp.y_);
56 Mat2::new(
57 Vector2::<f64>::new(s, -p),
58 Vector2::<f64>::new(-p * q, p * r + s),
59 )
60}
61
62/// The function `make_inverse` calculates the inverse of a 2x2 matrix.
63///
64/// Arguments:
65///
66/// * `vr`: A vector representing the row of a 2x2 matrix. The components of the vector are vr.x_ and
67/// vr.y_.
68/// * `vp`: The parameter `vp` represents a 2D vector with components `x` and `y`.
69///
70/// Returns:
71///
72/// The function `make_inverse` returns a `Mat2` object.
73#[inline]
74pub fn make_inverse(vr: &Vec2, vp: &Vec2) -> Mat2 {
75 let (r, q) = (vr.x_, vr.y_);
76 let (p, s) = (vp.x_, vp.y_);
77 let m_adjoint = Mat2::new(
78 Vector2::<f64>::new(s, -p),
79 Vector2::<f64>::new(-p * q, p * r + s),
80 );
81 m_adjoint / m_adjoint.det()
82}
83
84/// The `delta` function calculates the delta value for the Bairstow's method
85///
86/// Arguments:
87///
88/// * `vA`: A vector representing the coefficients of a polynomial equation.
89/// * `vr`: The parameter `vr` represents the vector `[-2.0, 0.0]`.
90/// * `vp`: The parameter `vp` represents the vector vr - vrj
91///
92/// Returns:
93///
94/// The function `delta` returns a `Vec2` object.
95///
96/// r * p - m -p
97/// q * p -m
98///
99/// # Examples:
100///
101/// ```
102/// use ginger::rootfinding::delta;
103/// use ginger::vector2::Vector2;
104///
105/// let mut vA1 = Vector2::new(1.0, 2.0);
106/// let vri = Vector2::new(-2.0, 0.0);
107/// let vrj = Vector2::new(4.0, 5.0);
108/// let vd = delta(&vA1, &vri, &vrj);
109/// assert_eq!(vd, Vector2::new(0.2, 0.4));
110/// ```
111#[inline]
112pub fn delta(vA: &Vec2, vr: &Vec2, vp: &Vec2) -> Vec2 {
113 let mp = make_adjoint(vr, vp); // 2 mul's
114 mp.mdot(vA) / mp.det() // 6 mul's + 2 div's
115}
116
117/// delta 1 for ri - rj
118///
119/// # Examples:
120///
121/// ```
122/// use ginger::rootfinding::delta1;
123/// use ginger::vector2::Vector2;
124///
125/// let mut vA1 = Vector2::new(1.0, 2.0);
126/// let vri = Vector2::new(-2.0, -0.0);
127/// let vrj = Vector2::new(4.0, -5.0);
128/// let vd = delta1(&vA1, &vri, &vrj);
129/// assert_eq!(vd, Vector2::new(0.2, 0.4));
130/// ```
131#[inline]
132pub fn delta1(vA: &Vec2, vr: &Vec2, vp: &Vec2) -> Vec2 {
133 let (r, q) = (vr.x_, vr.y_);
134 let (p, s) = (vp.x_, vp.y_);
135 let mp = Matrix2::new(Vec2::new(-s, -p), Vec2::new(p * q, p * r - s));
136 mp.mdot(vA) / mp.det() // 6 mul's + 2 div's
137}
138
139/// The `suppress_old` function performs zero suppression on a set of vectors.
140///
141/// Arguments:
142///
143/// * `vA`: A mutable reference to a Vector2 object representing the coefficients of a polynomial. The
144/// coefficients are stored in the x_ and y_ fields of the Vector2 object.
145/// * `vA1`: vA1 is a mutable reference to a Vector2 object.
146/// * `vri`: The parameter `vri` represents a vector with components `r` and `i`. It is used in the
147/// `suppress_old` function to perform calculations.
148/// * `vrj`: The parameter `vrj` represents a vector with components `x` and `y`.
149///
150/// # Examples:
151///
152/// ```
153/// use ginger::rootfinding::delta;
154/// use ginger::rootfinding::suppress_old;
155/// use ginger::vector2::Vector2;
156/// use approx_eq::assert_approx_eq;
157///
158/// let mut vA = Vector2::new(3.0, 3.0);
159/// let mut vA1 = Vector2::new(1.0, 2.0);
160/// let vri = Vector2::new(-2.0, 0.0);
161/// let vrj = Vector2::new(4.0, 5.0);
162///
163/// suppress_old(&mut vA, &mut vA1, &vri, &vrj);
164/// let dr = delta(&vA, &vri, &vA1);
165/// assert_approx_eq!(dr.x_, -16.780821917808325);
166/// assert_approx_eq!(dr.y_, 1.4383561643835612);
167#[inline]
168pub fn suppress_old(vA: &mut Vec2, vA1: &mut Vec2, vri: &Vec2, vrj: &Vec2) {
169 let (A, B) = (vA.x_, vA.y_);
170 let (A1, B1) = (vA1.x_, vA1.y_);
171 let vp = vri - vrj;
172 let (r, q) = (vri.x_, vri.y_);
173 let (p, s) = (vp.x_, vp.y_);
174 let f = (r * p) + s;
175 let qp = q * p;
176 let e = (f * s) - (qp * p);
177 let a = ((A * s) - (B * p)) / e;
178 let b = ((B * f) - (A * qp)) / e;
179 let c = A1 - a;
180 let d = (B1 - b) - (a * p);
181 vA.x_ = a;
182 vA.y_ = b;
183 vA1.x_ = ((c * s) - (d * p)) / e;
184 vA1.y_ = ((d * f) - (c * qp)) / e;
185}
186
187/// The `suppress` function in Rust performs zero suppression on a set of vectors.
188///
189/// Arguments:
190///
191/// * `vA`: A vector representing the coefficients of a polynomial function.
192/// * `vA1`: The parameter `vA1` is a `Vector2` object representing a vector with two components. It is
193/// used as an input parameter in the `suppress` function.
194/// * `vri`: The parameter `vri` represents the vector `ri`, and `vrj` represents the vector `rj`. These
195/// vectors are used in the calculation of the suppression step in the Bairstow's method for root
196/// finding.
197/// * `vrj`: The parameter `vrj` represents a vector with coordinates (4.0, 5.0).
198/// Zero suppression
199///
200/// # Examples:
201///
202/// ```
203/// use ginger::rootfinding::delta;
204/// use ginger::rootfinding::suppress;
205/// use ginger::vector2::Vector2;
206/// use approx_eq::assert_approx_eq;
207///
208/// let mut vA = Vector2::new(3.0, 3.0);
209/// let mut vA1 = Vector2::new(1.0, 2.0);
210/// let vri = Vector2::new(-2.0, 0.0);
211/// let vrj = Vector2::new(4.0, 5.0);
212///
213/// (vA, vA1) = suppress(&mut vA, &mut vA1, &vri, &vrj);
214/// let dr = delta(&vA, &vri, &vA1);
215/// assert_approx_eq!(dr.x_, -16.780821917808325);
216/// assert_approx_eq!(dr.y_, 1.4383561643835612);
217#[inline]
218pub fn suppress(vA: &Vec2, vA1: &Vec2, vri: &Vec2, vrj: &Vec2) -> (Vec2, Vec2) {
219 let vp = vri - vrj;
220 let m_inverse = make_inverse(vri, &vp);
221 let va = m_inverse.mdot(vA);
222 let mut vc = vA1 - va;
223 vc.y_ -= va.x_ * vp.x_;
224 let va1 = m_inverse.mdot(&vc);
225 (va, va1)
226}
227
228/// The `horner_eval` function in Rust implements the Horner's method for polynomial evaluation.
229///
230/// Arguments:
231///
232/// * `coeffs`: A mutable slice of f64 values representing the coefficients of a polynomial. The
233/// coefficients are ordered from highest degree to lowest degree.
234/// * `degree`: The `degree` parameter represents the degree of the polynomial. In the given example,
235/// the polynomial has a degree of 8.
236/// * `zval`: The `zval` parameter in the `horner_eval` function represents the value at which the
237/// polynomial is evaluated. It is the value of the independent variable in the polynomial expression.
238///
239/// Returns:
240///
241/// The function `horner_eval` returns a `f64` value, which is the result of evaluating the polynomial
242/// with the given coefficients at the specified value `zval`.
243///
244/// # Examples:
245///
246/// ```
247/// use ginger::rootfinding::horner_eval;
248/// use approx_eq::assert_approx_eq;
249///
250/// let coeffs = vec![10.0, 34.0, 75.0, 94.0, 150.0, 94.0, 75.0, 34.0, 10.0];
251/// let px = horner_eval(&coeffs, 2.0);
252///
253/// assert_approx_eq!(px, 18250.0);
254/// assert_approx_eq!(coeffs[3], 94.0);
255/// ```
256#[inline]
257pub fn horner_eval(coeffs: &[f64], zval: f64) -> f64 {
258 coeffs.iter().fold(0.0, |acc, coeff| acc * zval + coeff)
259}
260
261/// The `horner` function implements Horner's evaluation for Bairstow's method in Rust.
262///
263/// Arguments:
264///
265/// * `coeffs`: A mutable slice of f64 values representing the coefficients of the polynomial. The
266/// coefficients are in descending order of degree.
267/// * `degree`: The `degree` parameter represents the degree of the polynomial. It is used to determine
268/// the number of coefficients in the `coeffs` array.
269/// * `vr`: The parameter `vr` is a `Vec2` struct that contains two values, `x_` and `y_`. In the
270/// example, `vr` is initialized with the values `-1.0` and `-2.0`.
271///
272/// Returns:
273///
274/// The function `horner` returns a `Vec2` struct, which contains two `f64` values representing the
275/// results of the Horner evaluation.
276///
277/// # Examples:
278///
279/// ```
280/// use ginger::rootfinding::horner;
281/// use ginger::vector2::Vector2;
282/// use approx_eq::assert_approx_eq;
283///
284/// let mut coeffs = vec![10.0, 34.0, 75.0, 94.0, 150.0, 94.0, 75.0, 34.0, 10.0];
285/// let px = horner(&mut coeffs, 8, &Vector2::new(-1.0, -2.0));
286///
287/// assert_approx_eq!(px.x_, 114.0);
288/// assert_approx_eq!(px.y_, 134.0);
289/// assert_approx_eq!(coeffs[3], 15.0);
290/// ```
291pub fn horner(coeffs: &mut [f64], degree: usize, vr: &Vec2) -> Vec2 {
292 let Vec2 { x_: r, y_: q } = vr;
293 for idx in 0..(degree - 1) {
294 coeffs[idx + 1] += coeffs[idx] * r;
295 coeffs[idx + 2] += coeffs[idx] * q;
296 }
297 Vector2::<f64>::new(coeffs[degree - 1], coeffs[degree])
298}
299
300/// The `initial_guess` function in Rust calculates the initial guesses for the roots of a polynomial
301/// using Bairstow's method.
302///
303/// Arguments:
304///
305/// * `coeffs`: A vector of coefficients representing a polynomial.
306///
307/// Returns:
308///
309/// The function `initial_guess` returns a vector of `Vector2` structs, which represent the initial
310/// guesses for the roots of a polynomial equation.
311///
312/// # Examples:
313///
314/// ```
315/// use ginger::rootfinding::initial_guess;
316/// use ginger::vector2::Vector2;
317///
318/// let coeffs = vec![10.0, 34.0, 75.0, 94.0, 150.0, 94.0, 75.0, 34.0, 10.0];
319/// let vr0s = initial_guess(&coeffs);
320/// ```
321pub fn initial_guess(coeffs: &[f64]) -> Vec<Vec2> {
322 let mut degree = coeffs.len() - 1;
323 let center = -coeffs[1] / (coeffs[0] * degree as f64);
324 // let mut coeffs1 = coeffs.to_owned();
325 let centroid = horner_eval(coeffs, center); // ???
326 let re = centroid.abs().powf(1.0 / (degree as f64));
327 degree /= 2;
328 degree *= 2; // make even
329 let k = PI / (degree as f64);
330 let m = center * center + re * re;
331 (1..degree)
332 .step_by(2)
333 .map(|i| {
334 let temp = re * (k * i as f64).cos();
335 let r0 = 2.0 * (center + temp);
336 let t0 = m + 2.0 * center * temp;
337 Vector2::<f64>::new(r0, -t0)
338 })
339 .collect()
340}
341
342/// Parallel Bairstow's method (even degree only)
343///
344/// The `pbairstow_even` function implements the parallel Bairstow's method for finding roots of
345/// even-degree polynomials.
346///
347/// Arguments:
348///
349/// * `coeffs`: The `coeffs` parameter is a slice of `f64` values representing the coefficients of a polynomial.
350/// It is assumed that the polynomial has an even degree.
351/// * `vrs`: A vector of initial guesses for the roots of the polynomial. Each element of the vector is
352/// a complex number representing a root guess.
353/// * `options`: The `options` parameter is an instance of the `Options` struct, which contains the
354/// following fields:
355///
356/// # Examples:
357///
358/// ```
359/// use ginger::rootfinding::{initial_guess, pbairstow_even, Options};
360///
361/// let coeffs = vec![10.0, 34.0, 75.0, 94.0, 150.0, 94.0, 75.0, 34.0, 10.0];
362/// let mut vrs = initial_guess(&coeffs);
363/// let (niter, _found) = pbairstow_even(&coeffs, &mut vrs, &Options::default());
364///
365/// assert_eq!(niter, 5);
366/// ```
367pub fn pbairstow_even(coeffs: &[f64], vrs: &mut [Vec2], options: &Options) -> (usize, bool) {
368 let m_rs = vrs.len();
369 let mut converged = vec![false; m_rs];
370
371 for niter in 1..options.max_iters {
372 let mut tolerance = 0.0;
373 for i in 0..m_rs {
374 if converged[i] {
375 continue;
376 }
377 let mut vri = vrs[i];
378 if let Some(tol_i) = pbairstow_even_job(coeffs, i, &mut vri, &mut converged[i], vrs) {
379 if tolerance < tol_i {
380 tolerance = tol_i;
381 }
382 }
383 vrs[i] = vri;
384 }
385 if tolerance < options.tolerance {
386 return (niter, true);
387 }
388 }
389 (options.max_iters, false)
390}
391
392/// Multi-threading Bairstow's method (even degree only)
393///
394/// The `pbairstow_even_mt` function implements the multi-threading parallel Bairstow's
395/// method for finding roots of even-degree polynomials.
396///
397/// Arguments:
398///
399/// * `coeffs`: The `coeffs` parameter is a slice of `f64` values representing the coefficients of a polynomial.
400/// It is assumed that the polynomial has an even degree.
401/// * `vrs`: A vector of initial guesses for the roots of the polynomial. Each element of the vector is
402/// a complex number representing a root guess.
403/// * `options`: The `options` parameter is an instance of the `Options` struct, which contains the
404/// following fields:
405///
406/// # Examples:
407///
408/// ```
409/// use ginger::rootfinding::{initial_guess, pbairstow_even_mt, Options};
410///
411/// let coeffs = vec![10.0, 34.0, 75.0, 94.0, 150.0, 94.0, 75.0, 34.0, 10.0];
412/// let mut vrs = initial_guess(&coeffs);
413/// let (niter, _found) = pbairstow_even_mt(&coeffs, &mut vrs, &Options::default());
414///
415/// assert_eq!(niter, 8);
416/// ```
417pub fn pbairstow_even_mt(coeffs: &[f64], vrs: &mut Vec<Vec2>, options: &Options) -> (usize, bool) {
418 use rayon::prelude::*;
419
420 let m_rs = vrs.len();
421 let mut vrsc = vec![Vec2::default(); m_rs];
422 let mut converged = vec![false; m_rs];
423
424 for niter in 1..options.max_iters {
425 let mut tolerance = 0.0;
426 vrsc.copy_from_slice(vrs);
427
428 let tol_i = vrs
429 .par_iter_mut()
430 .zip(converged.par_iter_mut())
431 .enumerate()
432 .filter(|(_, (_, converged))| !**converged)
433 .filter_map(|(i, (vri, converged))| {
434 pbairstow_even_job(coeffs, i, vri, converged, &vrsc)
435 })
436 .reduce(|| tolerance, |x, y| x.max(y));
437 if tolerance < tol_i {
438 tolerance = tol_i;
439 }
440 if tolerance < options.tolerance {
441 return (niter, true);
442 }
443 }
444 (options.max_iters, false)
445}
446
447fn pbairstow_even_job(
448 coeffs: &[f64],
449 i: usize,
450 vri: &mut Vec2,
451 converged: &mut bool,
452 vrsc: &[Vec2],
453) -> Option<f64> {
454 let mut coeffs1 = coeffs.to_owned();
455 let degree = coeffs1.len() - 1; // degree, assume even
456 let mut vA = horner(&mut coeffs1, degree, vri);
457 let tol_i = vA.norm_inf();
458 if tol_i < 1e-15 {
459 *converged = true;
460 return None;
461 }
462 let mut vA1 = horner(&mut coeffs1, degree - 2, vri);
463 for (_, vrj) in vrsc.iter().enumerate().filter(|t| t.0 != i) {
464 suppress_old(&mut vA, &mut vA1, vri, vrj);
465 }
466 let dt = delta(&vA, vri, &vA1); // Gauss-Seidel fashion
467 *vri -= dt;
468 Some(tol_i)
469}
470
471/// The `initial_autocorr` function calculates the initial guesses for Bairstow's method for finding
472/// roots of a polynomial, specifically for the auto-correlation function.
473///
474/// Arguments:
475///
476/// * `coeffs`: The `coeffs` parameter is a slice of `f64` values representing the coefficients of a
477/// polynomial. The coefficients are ordered from highest degree to lowest degree.
478///
479/// Returns:
480///
481/// The function `initial_autocorr` returns a vector of `Vec2` structs.
482///
483/// # Examples:
484///
485/// ```
486/// use ginger::rootfinding::initial_autocorr;
487/// use ginger::vector2::Vector2;
488///
489/// let coeffs = vec![10.0, 34.0, 75.0, 94.0, 150.0, 94.0, 75.0, 34.0, 10.0];
490/// let vr0s = initial_autocorr(&coeffs);
491/// ```
492pub fn initial_autocorr(coeffs: &[f64]) -> Vec<Vec2> {
493 let mut degree = coeffs.len() - 1;
494 let re = coeffs[degree].abs().powf(1.0 / (degree as f64));
495 degree /= 2;
496 let k = PI / (degree as f64);
497 let m = re * re;
498 (1..degree)
499 .step_by(2)
500 .map(|i| Vector2::<f64>::new(2.0 * re * (k * i as f64).cos(), -m))
501 .collect()
502}
503
504/// The `pbairstow_autocorr` function implements the simultaneous Bairstow's method for finding roots of
505/// a polynomial, specifically for the auto-correlation function.
506///
507/// Arguments:
508///
509/// * `coeffs`: The `coeffs` parameter is a slice of `f64` values representing the coefficients of a
510/// polynomial. These coefficients are used to calculate the auto-correlation function.
511/// * `vrs`: `vrs` is a vector of complex numbers representing the initial guesses for the roots of the
512/// polynomial. Each element of `vrs` is a `Vec2` struct, which contains two fields: `x_` and `y_`.
513/// These fields represent the real and imaginary parts of the
514/// * `options`: The `Options` struct is used to specify the parameters for the Bairstow's method
515/// algorithm. It has the following fields:
516///
517/// # Examples:
518///
519/// ```
520/// use ginger::rootfinding::{initial_autocorr, pbairstow_autocorr, Options};
521///
522/// let coeffs = vec![10.0, 34.0, 75.0, 94.0, 150.0, 94.0, 75.0, 34.0, 10.0];
523/// let mut vrs = initial_autocorr(&coeffs);
524/// let (niter, _found) = pbairstow_autocorr(&coeffs, &mut vrs, &Options::default());
525///
526/// assert_eq!(niter, 1);
527/// ```
528pub fn pbairstow_autocorr(coeffs: &[f64], vrs: &mut [Vec2], options: &Options) -> (usize, bool) {
529 let m_rs = vrs.len();
530 let mut converged = vec![false; m_rs];
531
532 for niter in 0..options.max_iters {
533 let mut tolerance = 0.0;
534
535 for i in 0..m_rs {
536 if converged[i] {
537 continue;
538 }
539 let mut vri = vrs[i];
540 let tol_i = pbairstow_autocorr_mt_job(coeffs, i, &mut vri, &mut converged[i], vrs);
541 if let Some(tol_i) = tol_i {
542 if tolerance < tol_i {
543 tolerance = tol_i;
544 }
545 }
546 vrs[i] = vri;
547 }
548 if tolerance < options.tolerance {
549 return (niter, true);
550 }
551 }
552 (options.max_iters, false)
553}
554
555/// The `pbairstow_autocorr_mt` function is a multi-threaded implementation of Bairstow's method for
556/// finding roots of a polynomial, specifically for auto-correlation functions.
557///
558/// Arguments:
559///
560/// * `coeffs`: The `coeffs` parameter is a slice of `f64` values representing the coefficients of a
561/// polynomial. These coefficients are used as input for the Bairstow's method algorithm.
562/// * `vrs`: `vrs` is a vector of complex numbers representing the initial guesses for the roots of the
563/// polynomial. Each element of `vrs` is a `Vec2` struct, which contains the real and imaginary parts of
564/// the complex number.
565/// * `options`: The `options` parameter is an instance of the `Options` struct, which contains the
566/// following fields:
567///
568/// # Examples:
569///
570/// ```
571/// use ginger::rootfinding::{initial_autocorr, pbairstow_autocorr_mt, Options};
572///
573/// let coeffs = vec![10.0, 34.0, 75.0, 94.0, 150.0, 94.0, 75.0, 34.0, 10.0];
574/// let mut vrs = initial_autocorr(&coeffs);
575/// let (niter, _found) = pbairstow_autocorr_mt(&coeffs, &mut vrs, &Options::default());
576///
577/// assert_eq!(niter, 2);
578/// ```
579pub fn pbairstow_autocorr_mt(
580 coeffs: &[f64],
581 vrs: &mut Vec<Vec2>,
582 options: &Options,
583) -> (usize, bool) {
584 use rayon::prelude::*;
585
586 let m_rs = vrs.len();
587 let mut vrsc = vec![Vec2::default(); m_rs];
588 let mut converged = vec![false; m_rs];
589
590 for niter in 1..options.max_iters {
591 let mut tolerance = 0.0;
592 vrsc.copy_from_slice(vrs);
593
594 let tol_i = vrs
595 .par_iter_mut()
596 .zip(converged.par_iter_mut())
597 .enumerate()
598 .filter(|(_, (_, converged))| !**converged)
599 .filter_map(|(i, (vri, converged))| {
600 pbairstow_autocorr_mt_job(coeffs, i, vri, converged, &vrsc)
601 })
602 .reduce(|| tolerance, |x, y| x.max(y));
603 if tolerance < tol_i {
604 tolerance = tol_i;
605 }
606 if tolerance < options.tolerance {
607 return (niter, true);
608 }
609 }
610 (options.max_iters, false)
611}
612
613fn pbairstow_autocorr_mt_job(
614 coeffs: &[f64],
615 i: usize,
616 vri: &mut Vec2,
617 converged: &mut bool,
618 vrsc: &[Vec2],
619) -> Option<f64> {
620 let mut coeffs1 = coeffs.to_owned();
621 // let mut coeffs1 = coeffs.to_owned();
622 let degree = coeffs1.len() - 1; // assumed divided by 4
623 let mut vA = horner(&mut coeffs1, degree, vri);
624 let tol_i = vA.norm_inf();
625 if tol_i < 1e-15 {
626 *converged = true;
627 return None;
628 }
629 let mut vA1 = horner(&mut coeffs1, degree - 2, vri);
630 for (_j, vrj) in vrsc.iter().enumerate().filter(|t| t.0 != i) {
631 // vA1 -= delta(&vA, vrj, &(*vri - vrj));
632 suppress_old(&mut vA, &mut vA1, vri, vrj);
633 let vrjn = Vector2::<f64>::new(-vrj.x_, 1.0) / vrj.y_;
634 // vA1 -= delta(&vA, &vrjn, &(*vri - vrjn));
635 suppress_old(&mut vA, &mut vA1, vri, &vrjn);
636 }
637 let vrin = Vector2::<f64>::new(-vri.x_, 1.0) / vri.y_;
638 // vA1 -= delta(&vA, &vrin, &(*vri - vrin));
639 suppress_old(&mut vA, &mut vA1, vri, &vrin);
640 let dt = delta(&vA, vri, &vA1); // Gauss-Seidel fashion
641 *vri -= dt;
642 Some(tol_i)
643}
644
645/// The `extract_autocorr` function extracts the quadratic function where its roots are within a unit
646/// circle.
647///
648/// x^2 - r*x - t or x^2 + (r/t) * x + (-1/t)
649/// (x - a1)(x - a2) = x^2 - (a1 + a2) x + a1 * a2
650///
651/// Arguments:
652///
653/// * `vr`: A vector containing two values, representing the coefficients of a quadratic function. The
654/// first value represents the coefficient of x^2, and the second value represents the coefficient of x.
655///
656/// Returns:
657///
658/// The function `extract_autocorr` returns a `Vec2` struct, which contains two elements `x_` and `y_`.
659///
660/// # Examples:
661///
662/// ```
663/// use ginger::rootfinding::extract_autocorr;
664/// use ginger::vector2::Vector2;
665/// use approx_eq::assert_approx_eq;
666///
667/// let vr = extract_autocorr(Vector2::new(1.0, -4.0));
668///
669/// assert_approx_eq!(vr.x_, 0.25);
670/// assert_approx_eq!(vr.y_, -0.25);
671/// ```
672#[allow(dead_code)]
673pub fn extract_autocorr(vr: Vec2) -> Vec2 {
674 let Vec2 { x_: r, y_: q } = vr;
675 let hr = r / 2.0;
676 let d = hr * hr + q;
677 if d < 0.0 {
678 // complex conjugate root
679 if q < -1.0 {
680 return Vector2::<f64>::new(-r, 1.0) / q;
681 }
682 }
683 // two real roots
684 let mut a1 = hr + (if hr >= 0.0 { d.sqrt() } else { -d.sqrt() });
685 let mut a2 = -q / a1;
686
687 if a1.abs() > 1.0 {
688 if a2.abs() > 1.0 {
689 a2 = 1.0 / a2;
690 }
691 a1 = 1.0 / a1;
692 return Vector2::<f64>::new(a1 + a2, -a1 * a2);
693 }
694 if a2.abs() > 1.0 {
695 a2 = 1.0 / a2;
696 return Vector2::<f64>::new(a1 + a2, -a1 * a2);
697 }
698 // else no need to change
699 vr
700}