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
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
mod d;
mod non_central;
mod p;
mod q;
mod r;

use strafe_type::{LogProbability64, Positive64, Probability64, Rational64, Real64};

pub(crate) use self::{d::*, non_central::*, p::*, q::*, r::*};
use crate::traits::{Distribution, RNG};

/// # The F Distribution
///
/// ## Description:
///
/// Density, distribution function, quantile function and random
/// generation for the F distribution with ‘df1’ and ‘df2’ degrees of
/// freedom (and optional non-centrality parameter ‘ncp’).
///
/// ## Arguments:
///
/// * df1, df2: degrees of freedom.  ‘Inf’ is allowed.
/// * ncp: non-centrality parameter. If omitted the central F is
/// assumed.
///
/// ## Details:
///
/// The F distribution with ‘df1 =’ n1 and ‘df2 =’ n2 degrees of
/// freedom has density
///
/// $f(x) = \frac{\Gamma(\frac{n1 + n2}{2})}{\Gamma(\frac{n1}{2}) \Gamma(\frac{n2}{2})}
/// (\frac{n1}{n2})^{\frac{n1}{2}} x^{\frac{n1}{2} - 1}
/// (1 + \frac{n1}{n2} x)^{-\frac{n1 + n2}{2}}$
///
/// for $x > 0$.
///
/// It is the distribution of the ratio of the mean squares of n1 and
/// n2 independent standard normals, and hence of the ratio of two
/// independent chi-squared variates each divided by its degrees of
/// freedom.  Since the ratio of a normal and the root mean-square of
/// m independent normals has a Student's t_m distribution, the square
/// of a t_m variate has a F distribution on 1 and m degrees of
/// freedom.
///
/// The non-central F distribution is again the ratio of mean squares
/// of independent normals of unit variance, but those in the
/// numerator are allowed to have non-zero means and ‘ncp’ is the sum
/// of squares of the means.  See Chisquare for further details on
/// non-central distributions.
///
/// ## Value:
///
/// ‘df’ gives the density, ‘pf’ gives the distribution function ‘qf’
/// gives the quantile function, and ‘rf’ generates random deviates.
///
/// Invalid arguments will result in return value ‘NaN’, with a
/// warning.
///
/// The length of the result is determined by ‘n’ for ‘rf’, and is the
/// maximum of the lengths of the numerical arguments for the other
/// functions.
///
/// The numerical arguments other than ‘n’ are recycled to the length
/// of the result.  Only the first elements of the logical arguments
/// are used.
///
/// ## Density Plot
///
/// ```rust
/// # use r2rs_base::traits::StatisticalSlice;
/// # use r2rs_nmath::{distribution::FBuilder, traits::Distribution};
/// # use strafe_plot::prelude::{IntoDrawingArea, Line, Plot, PlotOptions, SVGBackend, BLACK};
/// # use strafe_type::FloatConstraint;
/// let f = FBuilder::new().build();
/// let x = <[f64]>::sequence(-0.5, 4.0, 1000);
/// let y = x
///     .iter()
///     .map(|x| f.density(x).unwrap())
///     .collect::<Vec<_>>();
///
/// let root = SVGBackend::new("density.svg", (1024, 768)).into_drawing_area();
/// Plot::new()
///     .with_options(PlotOptions {
///         x_axis_label: "x".to_string(),
///         y_axis_label: "density".to_string(),
///         ..Default::default()
///     })
///     .with_plottable(Line {
///         x,
///         y,
///         color: BLACK,
///         ..Default::default()
///     })
///     .plot(&root)
///     .unwrap();
/// # use std::fs::rename;
/// #     drop(root);
/// #     rename(
/// #             format!("density.svg"),
/// #             format!("src/distribution/f/doctest_out/density.svg"),
/// #     )
/// #     .unwrap();
/// ```
#[cfg_attr(feature = "doc_outputs", cfg_attr(all(), doc = embed_doc_image::embed_image!("density", "src/distribution/f/doctest_out/density.svg")))]
#[cfg_attr(feature = "doc_outputs", cfg_attr(all(), doc = "![Density][density]"))]
///
/// ## Note:
///
/// Supplying ‘ncp = 0’ uses the algorithm for the non-central
/// distribution, which is not the same algorithm used if ‘ncp’ is
/// omitted.  This is to give consistent behaviour in extreme cases
/// with values of ‘ncp’ very near zero.
///
/// The code for non-zero ‘ncp’ is principally intended to be used for
/// moderate values of ‘ncp’: it will not be highly accurate,
/// especially in the tails, for large values.
///
/// ## Source:
///
/// For the central case of ‘df’, computed _via_ a binomial
/// probability, code contributed by Catherine Loader (see ‘dbinom’);
/// for the non-central case computed _via_ ‘dbeta’, code contributed
/// by Peter Ruckdeschel.
///
/// For ‘pf’, _via_ ‘pbeta’ (or for large ‘df2’, _via_ ‘pchisq’).
///
/// For ‘qf’, _via_ ‘qchisq’ for large ‘df2’, else _via_ ‘qbeta’.
///
/// ## References:
///
/// Becker, R. A., Chambers, J. M. and Wilks, A. R. (1988) _The New S
/// Language_.  Wadsworth & Brooks/Cole.
///
/// Johnson, N. L., Kotz, S. and Balakrishnan, N. (1995) _Continuous
/// Univariate Distributions_, volume 2, chapters 27 and 30.  Wiley,
/// New York.
///
/// ## See Also:
///
/// Distributions for other standard distributions, including ‘dchisq’
/// for chi-squared and ‘dt’ for Student's t distributions.
///
/// ## Examples:
///
/// Lower Tails
/// ```rust
/// # use r2rs_base::traits::StatisticalSlice;
/// # use r2rs_nmath::{
/// #     distribution::{FBuilder, TBuilder},
/// #     traits::Distribution,
/// # };
/// # use strafe_type::FloatConstraint;
/// let x = <[f64]>::sequence(0.001, 5.0, 100);
/// let nu = 4.0;
///
/// let t = TBuilder::new().with_df(nu).unwrap().build();
/// let f = FBuilder::new().with_df1(1).with_df2(nu).build();
/// let r1 = x
///     .iter()
///     .map(|x| 2.0 * t.probability(x, true).unwrap() - 1.0)
///     .collect::<Vec<_>>();
/// let r2 = x
///     .iter()
///     .map(|x| f.probability(x.powi(2), true).unwrap())
///     .collect::<Vec<_>>();
/// println!("{r1:?}");
/// println!("{r2:?}");
/// # use std::{fs::File, io::Write};
/// # let mut f = File::create("src/distribution/f/doctest_out/lower_tail.md").unwrap();
/// # writeln!(f, "```output").unwrap();
/// # writeln!(f, "{r1:?}").unwrap();
/// # writeln!(f, "{r2:?}").unwrap();
/// # writeln!(f, "```").unwrap();
/// ```
#[cfg_attr(feature = "doc_outputs", cfg_attr(all(), doc = include_str!("doctest_out/lower_tail.md")))]
///
/// Upper Tails
/// ```rust
/// # use r2rs_base::traits::StatisticalSlice;
/// # use r2rs_nmath::{
/// #     distribution::{FBuilder, TBuilder},
/// #     traits::Distribution,
/// # };
/// # use strafe_type::FloatConstraint;
/// let x = <[f64]>::sequence(0.001, 5.0, 100);
/// let nu = 4.0;
///
/// let t = TBuilder::new().with_df(nu).unwrap().build();
/// let f = FBuilder::new().with_df1(1).with_df2(nu).build();
/// let r1 = x
///     .iter()
///     .map(|x| 2.0 * t.probability(x, false).unwrap())
///     .collect::<Vec<_>>();
/// let r2 = x
///     .iter()
///     .map(|x| f.probability(x.powi(2), false).unwrap())
///     .collect::<Vec<_>>();
/// println!("{r1:?}");
/// println!("{r2:?}");
/// # use std::{fs::File, io::Write};
/// # let mut f = File::create("src/distribution/f/doctest_out/upper_tail.md").unwrap();
/// # writeln!(f, "```output").unwrap();
/// # writeln!(f, "{r1:?}").unwrap();
/// # writeln!(f, "{r2:?}").unwrap();
/// # writeln!(f, "```").unwrap();
/// ```
#[cfg_attr(feature = "doc_outputs", cfg_attr(all(), doc = include_str!("doctest_out/upper_tail.md")))]
///
/// The density of the square of a $t_m$ is $2*\frac{\{t-density}(x, m)}{2*x}$
///
/// Check this is the same as the density of $F_{1,m}$
/// ```rust
/// # use r2rs_base::traits::StatisticalSlice;
/// # use r2rs_nmath::{
/// #     distribution::{FBuilder, TBuilder},
/// #     traits::Distribution,
/// # };
/// # use strafe_type::FloatConstraint;
/// let x = <[f64]>::sequence(0.001, 5.0, 100);
/// let nu = 5.0;
///
/// let t = TBuilder::new().with_df(nu).unwrap().build();
/// let f = FBuilder::new().with_df1(1).with_df2(nu).build();
/// let r1 = x
///     .iter()
///     .map(|x| t.density(x).unwrap() / x)
///     .collect::<Vec<_>>();
/// let r2 = x
///     .iter()
///     .map(|x| f.density(x.powi(2)).unwrap())
///     .collect::<Vec<_>>();
/// println!("{r1:?}");
/// println!("{r2:?}");
/// # use std::{fs::File, io::Write};
/// # let mut f = File::create("src/distribution/f/doctest_out/dens.md").unwrap();
/// # writeln!(f, "```output").unwrap();
/// # writeln!(f, "{r1:?}").unwrap();
/// # writeln!(f, "{r2:?}").unwrap();
/// # writeln!(f, "```").unwrap();
/// ```
#[cfg_attr(feature = "doc_outputs", cfg_attr(all(), doc = include_str!("doctest_out/dens.md")))]
///
/// Identity:  qf(2*p - 1, 1, df) == qt(p, df)^2  for  p >= 1/2
/// ```rust
/// use r2rs_base::traits::{QuantileType, StatisticalSlice};
/// use r2rs_nmath::{
///     distribution::{FBuilder, TBuilder},
///     traits::Distribution,
/// };
/// use strafe_type::FloatConstraint;
///
/// let p = <[f64]>::sequence(0.5, 0.99, 50);
/// let df = 10.0;
///
/// let rel_err = |x: &[f64], y: &[f64]| {
///     let mean = x
///         .iter()
///         .chain(y.iter())
///         .map(|f| f.abs())
///         .collect::<Vec<_>>()
///         .mean();
///     x.iter()
///         .zip(y.iter())
///         .map(|(x, y)| if x == y { 0.0 } else { (x - y).abs() / mean })
///         .collect::<Vec<_>>()
/// };
///
/// let t = TBuilder::new().with_df(df).unwrap().build();
/// let f = FBuilder::new().with_df1(1).with_df2(df).build();
/// let r1 = p
///     .iter()
///     .map(|p| t.quantile(p, true).unwrap().powi(2))
///     .collect::<Vec<_>>();
/// let r2 = p
///     .iter()
///     .map(|p| f.quantile(2.0 * p - 1.0, true).unwrap())
///     .collect::<Vec<_>>();
///
/// let error = rel_err(&r1, &r2).quantile(&[0.9], QuantileType::S);
///
/// println!("{error:?}");
/// # use std::{fs::File, io::Write};
/// # let mut f = File::create("src/distribution/f/doctest_out/quan.md").unwrap();
/// # writeln!(f, "```output").unwrap();
/// # writeln!(f, "{error:?}").unwrap();
/// # writeln!(f, "```").unwrap();
/// ```
#[cfg_attr(feature = "doc_outputs", cfg_attr(all(), doc = include_str!("doctest_out/quan.md")))]
pub struct F {
    df1: Rational64,
    df2: Rational64,
    ncp: Option<Positive64>,
}

impl Distribution for F {
    fn density<R: Into<Real64>>(&self, x: R) -> Real64 {
        if let Some(ncp) = self.ncp {
            dnf(x, self.df1, self.df2, ncp, false)
        } else {
            df(x, self.df1, self.df2, false)
        }
    }

    fn log_density<R: Into<Real64>>(&self, x: R) -> Real64 {
        if let Some(ncp) = self.ncp {
            dnf(x, self.df1, self.df2, ncp, true)
        } else {
            df(x, self.df1, self.df2, true)
        }
    }

    fn probability<R: Into<Real64>>(&self, q: R, lower_tail: bool) -> Probability64 {
        if let Some(ncp) = self.ncp {
            pnf(q, self.df1, self.df2, ncp, lower_tail)
        } else {
            pf(q, self.df1, self.df2, lower_tail)
        }
    }

    fn log_probability<R: Into<Real64>>(&self, q: R, lower_tail: bool) -> LogProbability64 {
        if let Some(ncp) = self.ncp {
            log_pnf(q, self.df1, self.df2, ncp, lower_tail)
        } else {
            log_pf(q, self.df1, self.df2, lower_tail)
        }
    }

    fn quantile<P: Into<Probability64>>(&self, p: P, lower_tail: bool) -> Real64 {
        if let Some(ncp) = self.ncp {
            qnf(p, self.df1, self.df2, ncp, lower_tail)
        } else {
            qf(p, self.df1, self.df2, lower_tail)
        }
    }

    fn log_quantile<LP: Into<LogProbability64>>(&self, p: LP, lower_tail: bool) -> Real64 {
        if let Some(ncp) = self.ncp {
            log_qnf(p, self.df1, self.df2, ncp, lower_tail)
        } else {
            log_qf(p, self.df1, self.df2, lower_tail)
        }
    }

    fn random_sample<R: RNG>(&self, rng: &mut R) -> Real64 {
        if let Some(ncp) = self.ncp {
            rnf(self.df1, self.df2, ncp, rng)
        } else {
            rf(self.df1, self.df2, rng)
        }
    }
}

pub struct FBuilder {
    df1: Option<Rational64>,
    df2: Option<Rational64>,
    ncp: Option<Positive64>,
}

impl FBuilder {
    pub fn new() -> Self {
        Self {
            df1: None,
            df2: None,
            ncp: None,
        }
    }

    pub fn with_df1<R: Into<Rational64>>(&mut self, df1: R) -> &mut Self {
        self.df1 = Some(df1.into());
        self
    }

    pub fn with_df2<R: Into<Rational64>>(&mut self, df2: R) -> &mut Self {
        self.df2 = Some(df2.into());
        self
    }

    pub fn with_ncp<P: Into<Positive64>>(&mut self, ncp: P) -> &mut Self {
        self.ncp = Some(ncp.into());
        self
    }

    pub fn build(&self) -> F {
        let df1 = self.df1.unwrap_or(1.0.into());
        let df2 = self.df2.unwrap_or(1.0.into());

        F {
            df1,
            df2,
            ncp: self.ncp,
        }
    }
}

#[cfg(test)]
mod tests;

#[cfg(all(test, feature = "enable_proptest"))]
mod proptests;

#[cfg(all(test, feature = "enable_covtest"))]
mod covtests;