Skip to main content

fdars_core/
scalar_on_function.rs

1//! Scalar-on-function regression with mixed scalar/functional covariates.
2//!
3//! Implements models of the form:
4//! ```text
5//! y = α + ∫β(t)X(t)dt + γᵀz + ε
6//! ```
7//! where X(t) is a functional predictor, z is a vector of scalar covariates,
8//! β(t) is the functional coefficient, and γ is the vector of scalar coefficients.
9//!
10//! # Methods
11//!
12//! - [`fregre_lm`]: FPC-based functional linear model with optional scalar covariates
13//! - [`fregre_np_mixed`]: Nonparametric kernel regression with product kernels
14//! - [`functional_logistic`]: Logistic regression for binary outcomes
15//! - [`fregre_cv`]: Cross-validation for number of FPC components
16
17use crate::cv::create_folds;
18use crate::iter_maybe_parallel;
19use crate::matrix::FdMatrix;
20use crate::regression::{fdata_to_pc_1d, FpcaResult};
21use rand::prelude::*;
22#[cfg(feature = "parallel")]
23use rayon::iter::ParallelIterator;
24
25// ---------------------------------------------------------------------------
26// Result types
27// ---------------------------------------------------------------------------
28
29/// Result of functional linear regression.
30pub struct FregreLmResult {
31    /// Intercept α
32    pub intercept: f64,
33    /// Functional coefficient β(t), evaluated on the original grid (length m)
34    pub beta_t: Vec<f64>,
35    /// Pointwise standard errors of β(t) (length m)
36    pub beta_se: Vec<f64>,
37    /// Scalar coefficients γ (one per scalar covariate)
38    pub gamma: Vec<f64>,
39    /// Fitted values ŷ (length n)
40    pub fitted_values: Vec<f64>,
41    /// Residuals y - ŷ (length n)
42    pub residuals: Vec<f64>,
43    /// R² statistic
44    pub r_squared: f64,
45    /// Adjusted R²
46    pub r_squared_adj: f64,
47    /// Standard errors of all coefficients (intercept, FPC scores, scalar covariates)
48    pub std_errors: Vec<f64>,
49    /// Number of FPC components used
50    pub ncomp: usize,
51    /// FPCA result (for projecting new data)
52    pub fpca: FpcaResult,
53    /// Regression coefficients on (FPC scores, scalar covariates) — internal
54    pub coefficients: Vec<f64>,
55    /// Residual standard error
56    pub residual_se: f64,
57    /// GCV criterion value (if computed)
58    pub gcv: f64,
59    /// Akaike Information Criterion
60    pub aic: f64,
61    /// Bayesian Information Criterion
62    pub bic: f64,
63}
64
65/// Result of nonparametric functional regression with mixed predictors.
66pub struct FregreNpResult {
67    /// Fitted values ŷ (length n)
68    pub fitted_values: Vec<f64>,
69    /// Residuals y - ŷ (length n)
70    pub residuals: Vec<f64>,
71    /// R² statistic
72    pub r_squared: f64,
73    /// Bandwidth for functional distance kernel
74    pub h_func: f64,
75    /// Bandwidth for scalar covariates kernel
76    pub h_scalar: f64,
77    /// Leave-one-out CV error
78    pub cv_error: f64,
79}
80
81/// Result of functional logistic regression.
82pub struct FunctionalLogisticResult {
83    /// Intercept α
84    pub intercept: f64,
85    /// Functional coefficient β(t), evaluated on the original grid (length m)
86    pub beta_t: Vec<f64>,
87    /// Pointwise standard errors of β(t) (length m)
88    pub beta_se: Vec<f64>,
89    /// Scalar coefficients γ (one per scalar covariate)
90    pub gamma: Vec<f64>,
91    /// Predicted probabilities P(Y=1) (length n)
92    pub probabilities: Vec<f64>,
93    /// Predicted class labels (0 or 1)
94    pub predicted_classes: Vec<u8>,
95    /// Number of FPC components used
96    pub ncomp: usize,
97    /// Classification accuracy on training data
98    pub accuracy: f64,
99    /// Standard errors of all coefficients (intercept, FPC scores, scalar covariates)
100    pub std_errors: Vec<f64>,
101    /// Regression coefficients on (FPC scores, scalar covariates) — internal
102    pub coefficients: Vec<f64>,
103    /// Log-likelihood at convergence
104    pub log_likelihood: f64,
105    /// Number of IRLS iterations
106    pub iterations: usize,
107    /// FPCA result (for projecting new data)
108    pub fpca: FpcaResult,
109    /// Akaike Information Criterion
110    pub aic: f64,
111    /// Bayesian Information Criterion
112    pub bic: f64,
113}
114
115/// Result of cross-validation for K selection.
116pub struct FregreCvResult {
117    /// Candidate K values tested
118    pub k_values: Vec<usize>,
119    /// CV error for each K
120    pub cv_errors: Vec<f64>,
121    /// Optimal K (minimizing CV error)
122    pub optimal_k: usize,
123    /// Minimum CV error
124    pub min_cv_error: f64,
125}
126
127/// Criterion used for model selection.
128#[derive(Debug, Clone, Copy, PartialEq)]
129pub enum SelectionCriterion {
130    /// Akaike Information Criterion
131    Aic,
132    /// Bayesian Information Criterion
133    Bic,
134    /// Generalized Cross-Validation
135    Gcv,
136}
137
138/// Result of ncomp model selection.
139pub struct ModelSelectionResult {
140    /// Best number of FPC components by the chosen criterion
141    pub best_ncomp: usize,
142    /// (ncomp, AIC, BIC, GCV) for each candidate
143    pub criteria: Vec<(usize, f64, f64, f64)>,
144}
145
146// ---------------------------------------------------------------------------
147// Shared linear algebra helpers
148// ---------------------------------------------------------------------------
149
150/// Compute X'X (symmetric, p×p stored flat row-major).
151pub(crate) fn compute_xtx(x: &FdMatrix) -> Vec<f64> {
152    let (n, p) = x.shape();
153    let mut xtx = vec![0.0; p * p];
154    for k in 0..p {
155        for j in k..p {
156            let mut s = 0.0;
157            for i in 0..n {
158                s += x[(i, k)] * x[(i, j)];
159            }
160            xtx[k * p + j] = s;
161            xtx[j * p + k] = s;
162        }
163    }
164    xtx
165}
166
167/// Compute X'y (length p).
168fn compute_xty(x: &FdMatrix, y: &[f64]) -> Vec<f64> {
169    let (n, p) = x.shape();
170    (0..p)
171        .map(|k| {
172            let mut s = 0.0;
173            for i in 0..n {
174                s += x[(i, k)] * y[i];
175            }
176            s
177        })
178        .collect()
179}
180
181/// Cholesky factorization: A = LL'. Returns L (p×p flat row-major) or None if singular.
182pub(crate) fn cholesky_factor(a: &[f64], p: usize) -> Option<Vec<f64>> {
183    let mut l = vec![0.0; p * p];
184    for j in 0..p {
185        let mut diag = a[j * p + j];
186        for k in 0..j {
187            diag -= l[j * p + k] * l[j * p + k];
188        }
189        if diag <= 1e-12 {
190            return None;
191        }
192        l[j * p + j] = diag.sqrt();
193        for i in (j + 1)..p {
194            let mut s = a[i * p + j];
195            for k in 0..j {
196                s -= l[i * p + k] * l[j * p + k];
197            }
198            l[i * p + j] = s / l[j * p + j];
199        }
200    }
201    Some(l)
202}
203
204/// Solve Lz = b (forward) then L'x = z (back). L is p×p flat row-major.
205pub(crate) fn cholesky_forward_back(l: &[f64], b: &[f64], p: usize) -> Vec<f64> {
206    let mut z = b.to_vec();
207    for j in 0..p {
208        for k in 0..j {
209            z[j] -= l[j * p + k] * z[k];
210        }
211        z[j] /= l[j * p + j];
212    }
213    for j in (0..p).rev() {
214        for k in (j + 1)..p {
215            z[j] -= l[k * p + j] * z[k];
216        }
217        z[j] /= l[j * p + j];
218    }
219    z
220}
221
222/// Solve Ax = b via Cholesky decomposition (A must be symmetric positive definite).
223fn cholesky_solve(a: &[f64], b: &[f64], p: usize) -> Option<Vec<f64>> {
224    let l = cholesky_factor(a, p)?;
225    Some(cholesky_forward_back(&l, b, p))
226}
227
228/// Compute hat matrix diagonal: H_ii = x_i' (X'X)^{-1} x_i, given Cholesky factor L of X'X.
229pub(crate) fn compute_hat_diagonal(x: &FdMatrix, l: &[f64]) -> Vec<f64> {
230    let (n, p) = x.shape();
231    let mut hat_diag = vec![0.0; n];
232    for i in 0..n {
233        let mut v = vec![0.0; p];
234        for j in 0..p {
235            v[j] = x[(i, j)];
236            for k in 0..j {
237                v[j] -= l[j * p + k] * v[k];
238            }
239            v[j] /= l[j * p + j];
240        }
241        hat_diag[i] = v.iter().map(|vi| vi * vi).sum();
242    }
243    hat_diag
244}
245
246/// Compute diagonal of (X'X)^{-1} given Cholesky factor L, then SE = sqrt(sigma² * diag).
247fn compute_ols_std_errors(l: &[f64], p: usize, sigma2: f64) -> Vec<f64> {
248    let mut se = vec![0.0; p];
249    for j in 0..p {
250        let mut v = vec![0.0; p];
251        v[j] = 1.0;
252        for k in 0..p {
253            for kk in 0..k {
254                v[k] -= l[k * p + kk] * v[kk];
255            }
256            v[k] /= l[k * p + k];
257        }
258        se[j] = (sigma2 * v.iter().map(|vi| vi * vi).sum::<f64>()).sqrt();
259    }
260    se
261}
262
263// ---------------------------------------------------------------------------
264// Design matrix and coefficient recovery
265// ---------------------------------------------------------------------------
266
267/// Build design matrix: \[1, ξ_1, ..., ξ_K, z_1, ..., z_p\].
268/// Validate inputs for fregre_lm / functional_logistic.
269fn validate_fregre_inputs(
270    n: usize,
271    m: usize,
272    y: &[f64],
273    scalar_covariates: Option<&FdMatrix>,
274) -> Option<()> {
275    if n < 3 || m == 0 || y.len() != n {
276        return None;
277    }
278    if let Some(sc) = scalar_covariates {
279        if sc.nrows() != n {
280            return None;
281        }
282    }
283    Some(())
284}
285
286/// Resolve ncomp: auto-select via CV if 0, otherwise clamp.
287fn resolve_ncomp(
288    ncomp: usize,
289    data: &FdMatrix,
290    y: &[f64],
291    scalar_covariates: Option<&FdMatrix>,
292    n: usize,
293    m: usize,
294) -> Option<usize> {
295    if ncomp == 0 {
296        let cv = fregre_cv(data, y, scalar_covariates, 1, m.min(n - 1).min(20), 5)?;
297        Some(cv.optimal_k)
298    } else {
299        Some(ncomp.min(n - 1).min(m))
300    }
301}
302
303pub(crate) fn build_design_matrix(
304    fpca_scores: &FdMatrix,
305    ncomp: usize,
306    scalar_covariates: Option<&FdMatrix>,
307    n: usize,
308) -> FdMatrix {
309    let p_scalar = scalar_covariates.map_or(0, |sc| sc.ncols());
310    let p_total = 1 + ncomp + p_scalar;
311    let mut design = FdMatrix::zeros(n, p_total);
312    for i in 0..n {
313        design[(i, 0)] = 1.0;
314        for k in 0..ncomp {
315            design[(i, 1 + k)] = fpca_scores[(i, k)];
316        }
317        if let Some(sc) = scalar_covariates {
318            for j in 0..p_scalar {
319                design[(i, 1 + ncomp + j)] = sc[(i, j)];
320            }
321        }
322    }
323    design
324}
325
326/// Recover functional coefficient β(t) = Σ_k γ_k φ_k(t).
327fn recover_beta_t(fpc_coeffs: &[f64], rotation: &FdMatrix, m: usize) -> Vec<f64> {
328    let ncomp = fpc_coeffs.len();
329    let mut beta_t = vec![0.0; m];
330    for k in 0..ncomp {
331        for j in 0..m {
332            beta_t[j] += fpc_coeffs[k] * rotation[(j, k)];
333        }
334    }
335    beta_t
336}
337
338/// Pointwise standard error of β(t) via error propagation through FPCA rotation.
339///
340/// SE[β(t_j)]² = Σ_k φ_k(t_j)² · SE[γ_k]²
341fn compute_beta_se(gamma_se: &[f64], rotation: &FdMatrix, m: usize) -> Vec<f64> {
342    let ncomp = gamma_se.len();
343    let mut beta_se = vec![0.0; m];
344    for j in 0..m {
345        let mut var_j = 0.0;
346        for k in 0..ncomp {
347            var_j += rotation[(j, k)].powi(2) * gamma_se[k].powi(2);
348        }
349        beta_se[j] = var_j.sqrt();
350    }
351    beta_se
352}
353
354/// Compute fitted values ŷ = X β.
355fn compute_fitted(design: &FdMatrix, coeffs: &[f64]) -> Vec<f64> {
356    let (n, p) = design.shape();
357    (0..n)
358        .map(|i| {
359            let mut yhat = 0.0;
360            for j in 0..p {
361                yhat += design[(i, j)] * coeffs[j];
362            }
363            yhat
364        })
365        .collect()
366}
367
368/// Compute R² and adjusted R².
369fn compute_r_squared(y: &[f64], residuals: &[f64], p_total: usize) -> (f64, f64) {
370    let n = y.len();
371    let y_mean = y.iter().sum::<f64>() / n as f64;
372    let ss_tot: f64 = y.iter().map(|&yi| (yi - y_mean).powi(2)).sum();
373    let ss_res: f64 = residuals.iter().map(|r| r * r).sum();
374    let r_squared = if ss_tot > 0.0 {
375        1.0 - ss_res / ss_tot
376    } else {
377        0.0
378    };
379    let df_model = (p_total - 1) as f64;
380    let r_squared_adj = if n as f64 - df_model - 1.0 > 0.0 {
381        1.0 - (1.0 - r_squared) * (n as f64 - 1.0) / (n as f64 - df_model - 1.0)
382    } else {
383        r_squared
384    };
385    (r_squared, r_squared_adj)
386}
387
388// ---------------------------------------------------------------------------
389// OLS solver
390// ---------------------------------------------------------------------------
391
392/// Solve ordinary least squares: min ||Xb - y||² via normal equations with Cholesky.
393/// Returns (coefficients, hat_matrix_diagonal) or None if singular.
394fn ols_solve(x: &FdMatrix, y: &[f64]) -> Option<(Vec<f64>, Vec<f64>)> {
395    let (n, p) = x.shape();
396    if n < p || p == 0 {
397        return None;
398    }
399    let xtx = compute_xtx(x);
400    let xty = compute_xty(x, y);
401    let l = cholesky_factor(&xtx, p)?;
402    let b = cholesky_forward_back(&l, &xty, p);
403    let hat_diag = compute_hat_diagonal(x, &l);
404    Some((b, hat_diag))
405}
406
407// ---------------------------------------------------------------------------
408// fregre_lm: FPC-based functional linear model
409// ---------------------------------------------------------------------------
410
411/// Functional linear model with optional scalar covariates.
412///
413/// Fits the model: `y = α + Σ_k γ_k ξ_k + γ_z' z + ε`
414/// where ξ_k are FPC scores of the functional predictor X(t) and z are scalar covariates.
415/// The functional coefficient is recovered as `β(t) = Σ_k γ_k φ_k(t)`.
416///
417/// # Arguments
418/// * `data` - Functional predictor matrix (n × m)
419/// * `y` - Scalar response vector (length n)
420/// * `scalar_covariates` - Optional scalar covariates matrix (n × p), `None` for pure functional model
421/// * `ncomp` - Number of FPC components (if 0, selected by GCV)
422///
423/// # Returns
424/// [`FregreLmResult`] with estimated coefficients, fitted values, and diagnostics
425pub fn fregre_lm(
426    data: &FdMatrix,
427    y: &[f64],
428    scalar_covariates: Option<&FdMatrix>,
429    ncomp: usize,
430) -> Option<FregreLmResult> {
431    let (n, m) = data.shape();
432    validate_fregre_inputs(n, m, y, scalar_covariates)?;
433
434    let ncomp = resolve_ncomp(ncomp, data, y, scalar_covariates, n, m)?;
435
436    let fpca = fdata_to_pc_1d(data, ncomp)?;
437    let design = build_design_matrix(&fpca.scores, ncomp, scalar_covariates, n);
438    let p_total = design.ncols();
439    let (coeffs, hat_diag) = ols_solve(&design, y)?;
440
441    let fitted_values = compute_fitted(&design, &coeffs);
442    let residuals: Vec<f64> = y
443        .iter()
444        .zip(&fitted_values)
445        .map(|(&yi, &yh)| yi - yh)
446        .collect();
447    let (r_squared, r_squared_adj) = compute_r_squared(y, &residuals, p_total);
448
449    let ss_res: f64 = residuals.iter().map(|r| r * r).sum();
450    let df_resid = (n as f64 - p_total as f64).max(1.0);
451    let residual_se = (ss_res / df_resid).sqrt();
452    let sigma2 = ss_res / df_resid;
453
454    let xtx = compute_xtx(&design);
455    let l = cholesky_factor(&xtx, p_total).unwrap_or_else(|| vec![1.0; p_total * p_total]);
456    let std_errors = compute_ols_std_errors(&l, p_total, sigma2);
457
458    let gcv = hat_diag
459        .iter()
460        .zip(&residuals)
461        .map(|(&h, &r)| (r / (1.0 - h).max(1e-10)).powi(2))
462        .sum::<f64>()
463        / n as f64;
464
465    let beta_t = recover_beta_t(&coeffs[1..1 + ncomp], &fpca.rotation, m);
466    let beta_se = compute_beta_se(&std_errors[1..1 + ncomp], &fpca.rotation, m);
467    let gamma: Vec<f64> = coeffs[1 + ncomp..].to_vec();
468
469    let nf = n as f64;
470    let rss = ss_res;
471    let aic = nf * (rss / nf).ln() + 2.0 * p_total as f64;
472    let bic = nf * (rss / nf).ln() + (nf).ln() * p_total as f64;
473
474    Some(FregreLmResult {
475        intercept: coeffs[0],
476        beta_t,
477        beta_se,
478        gamma,
479        fitted_values,
480        residuals,
481        r_squared,
482        r_squared_adj,
483        std_errors,
484        ncomp,
485        fpca,
486        coefficients: coeffs,
487        residual_se,
488        gcv,
489        aic,
490        bic,
491    })
492}
493
494// ---------------------------------------------------------------------------
495// fregre_cv: Cross-validation for K selection
496// ---------------------------------------------------------------------------
497
498/// Copy a subset of rows from src into dst.
499fn copy_rows(dst: &mut FdMatrix, src: &FdMatrix, src_rows: &[usize]) {
500    let ncols = src.ncols();
501    for (dst_i, &src_i) in src_rows.iter().enumerate() {
502        for j in 0..ncols {
503            dst[(dst_i, j)] = src[(src_i, j)];
504        }
505    }
506}
507
508/// Compute CV error for a single K across all folds (randomized assignment).
509fn cv_error_for_k(
510    data: &FdMatrix,
511    y: &[f64],
512    scalar_covariates: Option<&FdMatrix>,
513    k: usize,
514    n_folds: usize,
515    folds: &[usize],
516) -> Option<f64> {
517    let n = data.nrows();
518    let ncols = data.ncols();
519    let mut total_error = 0.0;
520    let mut count = 0;
521
522    for fold in 0..n_folds {
523        let train_idx: Vec<usize> = (0..n).filter(|&i| folds[i] != fold).collect();
524        let test_idx: Vec<usize> = (0..n).filter(|&i| folds[i] == fold).collect();
525        let n_test = test_idx.len();
526        if n_test == 0 || train_idx.len() < k + 2 {
527            continue;
528        }
529
530        let mut train_data = FdMatrix::zeros(train_idx.len(), ncols);
531        let mut test_data = FdMatrix::zeros(n_test, ncols);
532        copy_rows(&mut train_data, data, &train_idx);
533        copy_rows(&mut test_data, data, &test_idx);
534
535        let train_y: Vec<f64> = train_idx.iter().map(|&i| y[i]).collect();
536        let test_y: Vec<f64> = test_idx.iter().map(|&i| y[i]).collect();
537
538        let train_sc = scalar_covariates.map(|sc| {
539            let mut m = FdMatrix::zeros(train_idx.len(), sc.ncols());
540            copy_rows(&mut m, sc, &train_idx);
541            m
542        });
543        let test_sc = scalar_covariates.map(|sc| {
544            let mut m = FdMatrix::zeros(n_test, sc.ncols());
545            copy_rows(&mut m, sc, &test_idx);
546            m
547        });
548
549        let fit = match fregre_lm(&train_data, &train_y, train_sc.as_ref(), k) {
550            Some(f) => f,
551            None => continue,
552        };
553
554        let predictions = predict_fregre_lm(&fit, &test_data, test_sc.as_ref());
555        let fold_mse: f64 = predictions
556            .iter()
557            .zip(&test_y)
558            .map(|(&yhat, &yi)| (yhat - yi).powi(2))
559            .sum::<f64>()
560            / n_test as f64;
561
562        total_error += fold_mse * n_test as f64;
563        count += n_test;
564    }
565
566    if count > 0 {
567        Some(total_error / count as f64)
568    } else {
569        None
570    }
571}
572
573/// K-fold cross-validation for selecting the number of FPC components.
574///
575/// # Arguments
576/// * `data` - Functional predictor matrix (n × m)
577/// * `y` - Scalar response vector (length n)
578/// * `scalar_covariates` - Optional scalar covariates matrix
579/// * `k_min` - Minimum number of components to test
580/// * `k_max` - Maximum number of components to test
581/// * `n_folds` - Number of CV folds
582pub fn fregre_cv(
583    data: &FdMatrix,
584    y: &[f64],
585    scalar_covariates: Option<&FdMatrix>,
586    k_min: usize,
587    k_max: usize,
588    n_folds: usize,
589) -> Option<FregreCvResult> {
590    let n = data.nrows();
591    if n < n_folds || k_min < 1 || k_min > k_max {
592        return None;
593    }
594
595    let k_max = k_max.min(n - 2).min(data.ncols());
596
597    // Use randomized fold assignment (consistent seed for reproducibility)
598    let folds = create_folds(n, n_folds, 42);
599
600    let mut k_values = Vec::new();
601    let mut cv_errors = Vec::new();
602
603    for k in k_min..=k_max {
604        if let Some(err) = cv_error_for_k(data, y, scalar_covariates, k, n_folds, &folds) {
605            k_values.push(k);
606            cv_errors.push(err);
607        }
608    }
609
610    if k_values.is_empty() {
611        return None;
612    }
613
614    let (optimal_idx, &min_cv_error) = cv_errors
615        .iter()
616        .enumerate()
617        .min_by(|(_, a), (_, b)| a.partial_cmp(b).unwrap())
618        .unwrap();
619
620    Some(FregreCvResult {
621        k_values: k_values.clone(),
622        cv_errors,
623        optimal_k: k_values[optimal_idx],
624        min_cv_error,
625    })
626}
627
628// ---------------------------------------------------------------------------
629// Model selection
630// ---------------------------------------------------------------------------
631
632/// Select optimal ncomp for `fregre_lm` using AIC, BIC, or GCV.
633///
634/// Fits models for `ncomp = 1..=max_ncomp`, collects information criteria,
635/// and returns the best ncomp under the chosen criterion.
636///
637/// # Arguments
638/// * `data` — Functional predictor matrix (n × m)
639/// * `y` — Scalar responses (length n)
640/// * `scalar_covariates` — Optional scalar covariates (n × p)
641/// * `max_ncomp` — Maximum number of FPC components to try
642/// * `criterion` — Which criterion to minimise
643pub fn model_selection_ncomp(
644    data: &FdMatrix,
645    y: &[f64],
646    scalar_covariates: Option<&FdMatrix>,
647    max_ncomp: usize,
648    criterion: SelectionCriterion,
649) -> Option<ModelSelectionResult> {
650    if max_ncomp == 0 {
651        return None;
652    }
653
654    let mut criteria = Vec::with_capacity(max_ncomp);
655
656    for k in 1..=max_ncomp {
657        if let Some(fit) = fregre_lm(data, y, scalar_covariates, k) {
658            criteria.push((k, fit.aic, fit.bic, fit.gcv));
659        }
660    }
661
662    if criteria.is_empty() {
663        return None;
664    }
665
666    let best_idx = match criterion {
667        SelectionCriterion::Aic => criteria
668            .iter()
669            .enumerate()
670            .min_by(|(_, a), (_, b)| a.1.partial_cmp(&b.1).unwrap_or(std::cmp::Ordering::Equal))
671            .map(|(i, _)| i)
672            .unwrap_or(0),
673        SelectionCriterion::Bic => criteria
674            .iter()
675            .enumerate()
676            .min_by(|(_, a), (_, b)| a.2.partial_cmp(&b.2).unwrap_or(std::cmp::Ordering::Equal))
677            .map(|(i, _)| i)
678            .unwrap_or(0),
679        SelectionCriterion::Gcv => criteria
680            .iter()
681            .enumerate()
682            .min_by(|(_, a), (_, b)| a.3.partial_cmp(&b.3).unwrap_or(std::cmp::Ordering::Equal))
683            .map(|(i, _)| i)
684            .unwrap_or(0),
685    };
686
687    Some(ModelSelectionResult {
688        best_ncomp: criteria[best_idx].0,
689        criteria,
690    })
691}
692
693/// Predict new responses using a fitted functional linear model.
694///
695/// # Arguments
696/// * `fit` - A fitted [`FregreLmResult`]
697/// * `new_data` - New functional predictor matrix (n_new × m)
698/// * `new_scalar` - Optional new scalar covariates (n_new × p)
699pub fn predict_fregre_lm(
700    fit: &FregreLmResult,
701    new_data: &FdMatrix,
702    new_scalar: Option<&FdMatrix>,
703) -> Vec<f64> {
704    let (n_new, m) = new_data.shape();
705    let ncomp = fit.ncomp;
706    let p_scalar = fit.gamma.len();
707
708    let mut predictions = vec![0.0; n_new];
709    for i in 0..n_new {
710        let mut yhat = fit.intercept;
711        for k in 0..ncomp {
712            let mut s = 0.0;
713            for j in 0..m {
714                s += (new_data[(i, j)] - fit.fpca.mean[j]) * fit.fpca.rotation[(j, k)];
715            }
716            yhat += fit.coefficients[1 + k] * s;
717        }
718        if let Some(sc) = new_scalar {
719            for j in 0..p_scalar {
720                yhat += fit.gamma[j] * sc[(i, j)];
721            }
722        }
723        predictions[i] = yhat;
724    }
725    predictions
726}
727
728// ---------------------------------------------------------------------------
729// Nonparametric kernel regression
730// ---------------------------------------------------------------------------
731
732/// Gaussian kernel: K(d, h) = exp(-d² / (2h²))
733fn gaussian_kernel(d: f64, h: f64) -> f64 {
734    (-d * d / (2.0 * h * h)).exp()
735}
736
737/// Compute symmetric pairwise distance matrix (flat n×n).
738fn compute_pairwise_distances(data: &FdMatrix, argvals: &[f64]) -> Vec<f64> {
739    let n = data.nrows();
740    let weights = crate::helpers::simpsons_weights(argvals);
741    let mut dists = vec![0.0; n * n];
742    for i in 0..n {
743        for j in (i + 1)..n {
744            let d = crate::helpers::l2_distance(&data.row(i), &data.row(j), &weights);
745            dists[i * n + j] = d;
746            dists[j * n + i] = d;
747        }
748    }
749    dists
750}
751
752/// Compute pairwise Euclidean distance matrix for scalar covariates.
753fn compute_scalar_distances(sc: &FdMatrix) -> Vec<f64> {
754    let n = sc.nrows();
755    let p = sc.ncols();
756    let mut dists = vec![0.0; n * n];
757    for i in 0..n {
758        for j in (i + 1)..n {
759            let mut d2 = 0.0;
760            for k in 0..p {
761                let diff = sc[(i, k)] - sc[(j, k)];
762                d2 += diff * diff;
763            }
764            let d = d2.sqrt();
765            dists[i * n + j] = d;
766            dists[j * n + i] = d;
767        }
768    }
769    dists
770}
771
772/// Nadaraya-Watson LOO prediction for one observation.
773fn nw_loo_predict(
774    i: usize,
775    n: usize,
776    y: &[f64],
777    func_dists: &[f64],
778    scalar_dists: &[f64],
779    h_func: f64,
780    h_scalar: f64,
781    has_scalar: bool,
782) -> f64 {
783    let mut num = 0.0;
784    let mut den = 0.0;
785    for j in 0..n {
786        if i == j {
787            continue;
788        }
789        let kf = gaussian_kernel(func_dists[i * n + j], h_func);
790        let ks = if has_scalar {
791            gaussian_kernel(scalar_dists[i * n + j], h_scalar)
792        } else {
793            1.0
794        };
795        let w = kf * ks;
796        num += w * y[j];
797        den += w;
798    }
799    if den > 1e-15 {
800        num / den
801    } else {
802        y[i]
803    }
804}
805
806/// LOO-CV error for Nadaraya-Watson with a single bandwidth.
807fn loo_cv_error(dists: &[f64], y: &[f64], n: usize, h: f64) -> f64 {
808    (0..n)
809        .map(|i| {
810            let mut num = 0.0;
811            let mut den = 0.0;
812            for j in 0..n {
813                if i == j {
814                    continue;
815                }
816                let w = gaussian_kernel(dists[i * n + j], h);
817                num += w * y[j];
818                den += w;
819            }
820            let yhat = if den > 1e-15 { num / den } else { y[i] };
821            (y[i] - yhat).powi(2)
822        })
823        .sum::<f64>()
824        / n as f64
825}
826
827/// Select bandwidth by minimizing LOO-CV error on a grid of distance quantiles.
828fn select_bandwidth_loo(dists: &[f64], y: &[f64], n: usize, _other_dists: Option<&[f64]>) -> f64 {
829    let mut nonzero_dists: Vec<f64> = (0..n)
830        .flat_map(|i| ((i + 1)..n).map(move |j| dists[i * n + j]))
831        .filter(|&d| d > 0.0)
832        .collect();
833    nonzero_dists.sort_by(|a, b| a.partial_cmp(b).unwrap());
834
835    if nonzero_dists.is_empty() {
836        return 1.0;
837    }
838
839    let n_cand = 20;
840    let mut best_h = nonzero_dists[nonzero_dists.len() / 2];
841    let mut best_cv = f64::INFINITY;
842
843    for qi in 1..=n_cand {
844        let q = qi as f64 / (n_cand + 1) as f64;
845        let idx = ((nonzero_dists.len() as f64 * q) as usize).min(nonzero_dists.len() - 1);
846        let h = nonzero_dists[idx].max(1e-10);
847        let cv = loo_cv_error(dists, y, n, h);
848        if cv < best_cv {
849            best_cv = cv;
850            best_h = h;
851        }
852    }
853    best_h
854}
855
856/// Nonparametric kernel regression with mixed functional and scalar predictors.
857///
858/// Uses product kernels:
859/// ```text
860/// ŷ(x, z) = Σᵢ K_func(Xᵢ, x) · K_scalar(zᵢ, z) · yᵢ / Σᵢ K_func(Xᵢ, x) · K_scalar(zᵢ, z)
861/// ```
862///
863/// Bandwidths are selected via leave-one-out CV if set to 0.
864///
865/// # Arguments
866/// * `data` - Functional predictor matrix (n × m)
867/// * `y` - Scalar response vector
868/// * `argvals` - Grid points for integration (length m)
869/// * `scalar_covariates` - Optional scalar covariates (n × p)
870/// * `h_func` - Bandwidth for functional kernel (0 for automatic)
871/// * `h_scalar` - Bandwidth for scalar kernel (0 for automatic)
872pub fn fregre_np_mixed(
873    data: &FdMatrix,
874    y: &[f64],
875    argvals: &[f64],
876    scalar_covariates: Option<&FdMatrix>,
877    h_func: f64,
878    h_scalar: f64,
879) -> Option<FregreNpResult> {
880    let n = data.nrows();
881    if n < 3 || data.ncols() == 0 || y.len() != n || argvals.len() != data.ncols() {
882        return None;
883    }
884
885    let func_dists = compute_pairwise_distances(data, argvals);
886    let has_scalar = scalar_covariates.is_some();
887    let scalar_dists = scalar_covariates
888        .map(compute_scalar_distances)
889        .unwrap_or_default();
890
891    let h_func = if h_func <= 0.0 {
892        select_bandwidth_loo(&func_dists, y, n, None)
893    } else {
894        h_func
895    };
896
897    let h_scalar = if has_scalar && h_scalar <= 0.0 {
898        select_bandwidth_loo(&scalar_dists, y, n, Some(&func_dists))
899    } else {
900        h_scalar
901    };
902
903    let mut fitted_values = vec![0.0; n];
904    let mut cv_error = 0.0;
905    for i in 0..n {
906        fitted_values[i] = nw_loo_predict(
907            i,
908            n,
909            y,
910            &func_dists,
911            &scalar_dists,
912            h_func,
913            h_scalar,
914            has_scalar,
915        );
916        cv_error += (y[i] - fitted_values[i]).powi(2);
917    }
918    cv_error /= n as f64;
919
920    let residuals: Vec<f64> = y
921        .iter()
922        .zip(&fitted_values)
923        .map(|(&yi, &yh)| yi - yh)
924        .collect();
925    let (r_squared, _) = compute_r_squared(y, &residuals, 1);
926
927    Some(FregreNpResult {
928        fitted_values,
929        residuals,
930        r_squared,
931        h_func,
932        h_scalar,
933        cv_error,
934    })
935}
936
937/// Predict new responses using a fitted nonparametric model.
938pub fn predict_fregre_np(
939    train_data: &FdMatrix,
940    y: &[f64],
941    train_scalar: Option<&FdMatrix>,
942    new_data: &FdMatrix,
943    new_scalar: Option<&FdMatrix>,
944    argvals: &[f64],
945    h_func: f64,
946    h_scalar: f64,
947) -> Vec<f64> {
948    let n_train = train_data.nrows();
949    let n_new = new_data.nrows();
950    let weights = crate::helpers::simpsons_weights(argvals);
951
952    (0..n_new)
953        .map(|i| {
954            let new_row = new_data.row(i);
955            let mut num = 0.0;
956            let mut den = 0.0;
957            for j in 0..n_train {
958                let d_func = crate::helpers::l2_distance(&new_row, &train_data.row(j), &weights);
959                let kf = gaussian_kernel(d_func, h_func);
960                let ks = match (new_scalar, train_scalar) {
961                    (Some(ns), Some(ts)) => {
962                        let d2: f64 = (0..ns.ncols())
963                            .map(|k| (ns[(i, k)] - ts[(j, k)]).powi(2))
964                            .sum();
965                        gaussian_kernel(d2.sqrt(), h_scalar)
966                    }
967                    _ => 1.0,
968                };
969                let w = kf * ks;
970                num += w * y[j];
971                den += w;
972            }
973            if den > 1e-15 {
974                num / den
975            } else {
976                0.0
977            }
978        })
979        .collect()
980}
981
982// ---------------------------------------------------------------------------
983// Functional logistic regression
984// ---------------------------------------------------------------------------
985
986/// One IRLS step: compute working response and solve weighted least squares.
987/// Returns updated beta or None if system is singular.
988fn irls_step(design: &FdMatrix, y: &[f64], beta: &[f64]) -> Option<Vec<f64>> {
989    let (n, p) = design.shape();
990
991    // Linear predictor η = Xβ, probabilities μ = sigmoid(η)
992    let eta: Vec<f64> = (0..n)
993        .map(|i| (0..p).map(|j| design[(i, j)] * beta[j]).sum())
994        .collect();
995    let mu: Vec<f64> = eta.iter().map(|&e| sigmoid(e)).collect();
996    let w: Vec<f64> = mu.iter().map(|&p| (p * (1.0 - p)).max(1e-10)).collect();
997    let z_work: Vec<f64> = (0..n).map(|i| eta[i] + (y[i] - mu[i]) / w[i]).collect();
998
999    // Weighted normal equations: (X'WX) β = X'Wz
1000    let mut xtwx = vec![0.0; p * p];
1001    for k in 0..p {
1002        for j in k..p {
1003            let mut s = 0.0;
1004            for i in 0..n {
1005                s += design[(i, k)] * w[i] * design[(i, j)];
1006            }
1007            xtwx[k * p + j] = s;
1008            xtwx[j * p + k] = s;
1009        }
1010    }
1011
1012    let mut xtwz = vec![0.0; p];
1013    for k in 0..p {
1014        let mut s = 0.0;
1015        for i in 0..n {
1016            s += design[(i, k)] * w[i] * z_work[i];
1017        }
1018        xtwz[k] = s;
1019    }
1020
1021    cholesky_solve(&xtwx, &xtwz, p)
1022}
1023
1024/// Compute log-likelihood of logistic model.
1025fn logistic_log_likelihood(probabilities: &[f64], y: &[f64]) -> f64 {
1026    probabilities
1027        .iter()
1028        .zip(y)
1029        .map(|(&p, &yi)| {
1030            let p = p.clamp(1e-15, 1.0 - 1e-15);
1031            yi * p.ln() + (1.0 - yi) * (1.0 - p).ln()
1032        })
1033        .sum()
1034}
1035
1036/// Sigmoid function: 1 / (1 + exp(-x))
1037pub(crate) fn sigmoid(x: f64) -> f64 {
1038    if x >= 0.0 {
1039        1.0 / (1.0 + (-x).exp())
1040    } else {
1041        let ex = x.exp();
1042        ex / (1.0 + ex)
1043    }
1044}
1045
1046/// Run IRLS iteration loop and return (beta, iterations).
1047fn irls_loop(design: &FdMatrix, y: &[f64], max_iter: usize, tol: f64) -> (Vec<f64>, usize) {
1048    let p_total = design.ncols();
1049    let mut beta = vec![0.0; p_total];
1050    let mut iterations = 0;
1051    for iter in 0..max_iter {
1052        iterations = iter + 1;
1053        let beta_new = match irls_step(design, y, &beta) {
1054            Some(b) => b,
1055            None => break,
1056        };
1057        let change: f64 = beta_new
1058            .iter()
1059            .zip(&beta)
1060            .map(|(a, b)| (a - b).abs())
1061            .fold(0.0_f64, f64::max);
1062        beta = beta_new;
1063        if change < tol {
1064            break;
1065        }
1066    }
1067    (beta, iterations)
1068}
1069
1070/// Build logistic result from converged beta.
1071fn build_logistic_result(
1072    design: &FdMatrix,
1073    beta: Vec<f64>,
1074    y: &[f64],
1075    fpca: FpcaResult,
1076    ncomp: usize,
1077    m: usize,
1078    iterations: usize,
1079) -> FunctionalLogisticResult {
1080    let (n, p) = design.shape();
1081    let eta = compute_fitted(design, &beta);
1082    let probabilities: Vec<f64> = eta.iter().map(|&e| sigmoid(e)).collect();
1083    let predicted_classes: Vec<u8> = probabilities
1084        .iter()
1085        .map(|&p| if p >= 0.5 { 1 } else { 0 })
1086        .collect();
1087    let correct: usize = predicted_classes
1088        .iter()
1089        .zip(y)
1090        .filter(|(&pred, &actual)| pred as f64 == actual)
1091        .count();
1092    let beta_t = recover_beta_t(&beta[1..1 + ncomp], &fpca.rotation, m);
1093    let gamma: Vec<f64> = beta[1 + ncomp..].to_vec();
1094
1095    // Compute coefficient SEs from Fisher information matrix (X'WX)^{-1}
1096    let w: Vec<f64> = probabilities
1097        .iter()
1098        .map(|&mu| (mu * (1.0 - mu)).max(1e-10))
1099        .collect();
1100    let mut xtwx = vec![0.0; p * p];
1101    for k in 0..p {
1102        for j in k..p {
1103            let mut s = 0.0;
1104            for i in 0..n {
1105                s += design[(i, k)] * w[i] * design[(i, j)];
1106            }
1107            xtwx[k * p + j] = s;
1108            xtwx[j * p + k] = s;
1109        }
1110    }
1111    let std_errors = cholesky_factor(&xtwx, p)
1112        .map(|l| compute_ols_std_errors(&l, p, 1.0))
1113        .unwrap_or_else(|| vec![f64::NAN; p]);
1114    let beta_se = compute_beta_se(&std_errors[1..1 + ncomp], &fpca.rotation, m);
1115
1116    let ll = logistic_log_likelihood(&probabilities, y);
1117    let deviance = -2.0 * ll;
1118    let nf = n as f64;
1119    let pf = p as f64;
1120    let aic = deviance + 2.0 * pf;
1121    let bic = deviance + nf.ln() * pf;
1122
1123    FunctionalLogisticResult {
1124        intercept: beta[0],
1125        beta_t,
1126        beta_se,
1127        gamma,
1128        accuracy: correct as f64 / nf,
1129        log_likelihood: ll,
1130        probabilities,
1131        predicted_classes,
1132        ncomp,
1133        std_errors,
1134        coefficients: beta,
1135        iterations,
1136        fpca,
1137        aic,
1138        bic,
1139    }
1140}
1141
1142/// Functional logistic regression for binary outcomes.
1143///
1144/// Fits: `log(P(Y=1)/P(Y=0)) = α + ∫β(t)X(t)dt + γᵀz`
1145/// using IRLS (iteratively reweighted least squares) on FPC scores.
1146///
1147/// # Arguments
1148/// * `data` - Functional predictor matrix (n × m)
1149/// * `y` - Binary response vector (0.0 or 1.0, length n)
1150/// * `scalar_covariates` - Optional scalar covariates (n × p)
1151/// * `ncomp` - Number of FPC components
1152/// * `max_iter` - Maximum IRLS iterations (default: 25)
1153/// * `tol` - Convergence tolerance (default: 1e-6)
1154pub fn functional_logistic(
1155    data: &FdMatrix,
1156    y: &[f64],
1157    scalar_covariates: Option<&FdMatrix>,
1158    ncomp: usize,
1159    max_iter: usize,
1160    tol: f64,
1161) -> Option<FunctionalLogisticResult> {
1162    let (n, m) = data.shape();
1163    if n < 3 || m == 0 || y.len() != n {
1164        return None;
1165    }
1166    if y.iter().any(|&yi| yi != 0.0 && yi != 1.0) {
1167        return None;
1168    }
1169
1170    let ncomp = ncomp.min(n - 1).min(m);
1171    let fpca = fdata_to_pc_1d(data, ncomp)?;
1172    let design = build_design_matrix(&fpca.scores, ncomp, scalar_covariates, n);
1173
1174    let max_iter = if max_iter == 0 { 25 } else { max_iter };
1175    let tol = if tol <= 0.0 { 1e-6 } else { tol };
1176
1177    let (beta, iterations) = irls_loop(&design, y, max_iter, tol);
1178    Some(build_logistic_result(
1179        &design, beta, y, fpca, ncomp, m, iterations,
1180    ))
1181}
1182
1183// ---------------------------------------------------------------------------
1184// Bootstrap CIs for β(t)
1185// ---------------------------------------------------------------------------
1186
1187/// Result of bootstrap confidence intervals for β(t).
1188pub struct BootstrapCiResult {
1189    /// Pointwise lower bound (length m).
1190    pub lower: Vec<f64>,
1191    /// Pointwise upper bound (length m).
1192    pub upper: Vec<f64>,
1193    /// Original β(t) estimate (length m).
1194    pub center: Vec<f64>,
1195    /// Simultaneous lower bound (sup-norm adjusted, length m).
1196    pub sim_lower: Vec<f64>,
1197    /// Simultaneous upper bound (sup-norm adjusted, length m).
1198    pub sim_upper: Vec<f64>,
1199    /// Number of bootstrap replicates that converged.
1200    pub n_boot_success: usize,
1201}
1202
1203/// Gather rows from `src` by index (with replacement), returning a new matrix.
1204fn subsample_rows(src: &FdMatrix, indices: &[usize]) -> FdMatrix {
1205    let ncols = src.ncols();
1206    let mut out = FdMatrix::zeros(indices.len(), ncols);
1207    for (dst_i, &src_i) in indices.iter().enumerate() {
1208        for j in 0..ncols {
1209            out[(dst_i, j)] = src[(src_i, j)];
1210        }
1211    }
1212    out
1213}
1214
1215/// Compute the q-th quantile of a sorted slice (linear interpolation).
1216fn quantile(sorted: &[f64], q: f64) -> f64 {
1217    if sorted.is_empty() {
1218        return f64::NAN;
1219    }
1220    if sorted.len() == 1 {
1221        return sorted[0];
1222    }
1223    let pos = q * (sorted.len() - 1) as f64;
1224    let lo = pos.floor() as usize;
1225    let hi = lo + 1;
1226    let frac = pos - lo as f64;
1227    if hi >= sorted.len() {
1228        sorted[sorted.len() - 1]
1229    } else {
1230        sorted[lo] * (1.0 - frac) + sorted[hi] * frac
1231    }
1232}
1233
1234/// Bootstrap confidence intervals for β(t) from a functional linear model.
1235///
1236/// Uses cases bootstrap (resampling observation indices with replacement) to
1237/// build pointwise and simultaneous confidence bands for the functional
1238/// coefficient β(t).
1239///
1240/// # Arguments
1241/// * `data` — Functional predictor matrix (n × m)
1242/// * `y` — Scalar response vector (length n)
1243/// * `scalar_covariates` — Optional scalar covariates (n × p)
1244/// * `ncomp` — Number of FPC components
1245/// * `n_boot` — Number of bootstrap replicates
1246/// * `alpha` — Significance level (e.g., 0.05 for 95% CI)
1247/// * `seed` — RNG seed for reproducibility
1248pub fn bootstrap_ci_fregre_lm(
1249    data: &FdMatrix,
1250    y: &[f64],
1251    scalar_covariates: Option<&FdMatrix>,
1252    ncomp: usize,
1253    n_boot: usize,
1254    alpha: f64,
1255    seed: u64,
1256) -> Option<BootstrapCiResult> {
1257    let (n, m) = data.shape();
1258    if n < 3 || m == 0 || y.len() != n || n_boot == 0 || alpha <= 0.0 || alpha >= 1.0 {
1259        return None;
1260    }
1261
1262    // Fit original model
1263    let original = fregre_lm(data, y, scalar_covariates, ncomp)?;
1264    let center = original.beta_t.clone();
1265
1266    // Bootstrap replicates
1267    let boot_betas: Vec<Vec<f64>> = iter_maybe_parallel!(0..n_boot)
1268        .filter_map(|b| {
1269            let mut rng = StdRng::seed_from_u64(seed.wrapping_add(b as u64));
1270            let indices: Vec<usize> = (0..n).map(|_| rng.gen_range(0..n)).collect();
1271
1272            let boot_data = subsample_rows(data, &indices);
1273            let boot_y: Vec<f64> = indices.iter().map(|&i| y[i]).collect();
1274            let boot_sc = scalar_covariates.map(|sc| subsample_rows(sc, &indices));
1275
1276            fregre_lm(&boot_data, &boot_y, boot_sc.as_ref(), ncomp).map(|fit| fit.beta_t)
1277        })
1278        .collect();
1279
1280    let n_boot_success = boot_betas.len();
1281    if n_boot_success < 3 {
1282        return None;
1283    }
1284
1285    // Pointwise bands: sort each grid point across replicates
1286    let lo_q = alpha / 2.0;
1287    let hi_q = 1.0 - alpha / 2.0;
1288    let mut lower = vec![0.0; m];
1289    let mut upper = vec![0.0; m];
1290    let mut boot_se = vec![0.0; m];
1291
1292    for j in 0..m {
1293        let mut vals: Vec<f64> = boot_betas.iter().map(|b| b[j]).collect();
1294        vals.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
1295        lower[j] = quantile(&vals, lo_q);
1296        upper[j] = quantile(&vals, hi_q);
1297
1298        // Bootstrap SE at this grid point
1299        let mean_j: f64 = vals.iter().sum::<f64>() / vals.len() as f64;
1300        let var_j: f64 =
1301            vals.iter().map(|&v| (v - mean_j).powi(2)).sum::<f64>() / (vals.len() - 1) as f64;
1302        boot_se[j] = var_j.sqrt().max(1e-15);
1303    }
1304
1305    // Simultaneous bands: sup-norm bootstrap
1306    let mut t_stats: Vec<f64> = boot_betas
1307        .iter()
1308        .map(|b| {
1309            (0..m)
1310                .map(|j| ((b[j] - center[j]) / boot_se[j]).abs())
1311                .fold(0.0_f64, f64::max)
1312        })
1313        .collect();
1314    t_stats.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
1315
1316    let c_alpha = quantile(&t_stats, 1.0 - alpha);
1317    let sim_lower: Vec<f64> = (0..m).map(|j| center[j] - c_alpha * boot_se[j]).collect();
1318    let sim_upper: Vec<f64> = (0..m).map(|j| center[j] + c_alpha * boot_se[j]).collect();
1319
1320    Some(BootstrapCiResult {
1321        lower,
1322        upper,
1323        center,
1324        sim_lower,
1325        sim_upper,
1326        n_boot_success,
1327    })
1328}
1329
1330/// Bootstrap confidence intervals for β(t) from a functional logistic model.
1331///
1332/// Same algorithm as [`bootstrap_ci_fregre_lm`] but using [`functional_logistic`]
1333/// as the inner estimator. Degenerate resamples (all-0 or all-1 y) fail naturally.
1334///
1335/// # Arguments
1336/// * `data` — Functional predictor matrix (n × m)
1337/// * `y` — Binary response vector (0.0 or 1.0, length n)
1338/// * `scalar_covariates` — Optional scalar covariates (n × p)
1339/// * `ncomp` — Number of FPC components
1340/// * `n_boot` — Number of bootstrap replicates
1341/// * `alpha` — Significance level
1342/// * `seed` — RNG seed
1343/// * `max_iter` — Maximum IRLS iterations per replicate
1344/// * `tol` — IRLS convergence tolerance
1345pub fn bootstrap_ci_functional_logistic(
1346    data: &FdMatrix,
1347    y: &[f64],
1348    scalar_covariates: Option<&FdMatrix>,
1349    ncomp: usize,
1350    n_boot: usize,
1351    alpha: f64,
1352    seed: u64,
1353    max_iter: usize,
1354    tol: f64,
1355) -> Option<BootstrapCiResult> {
1356    let (n, m) = data.shape();
1357    if n < 3 || m == 0 || y.len() != n || n_boot == 0 || alpha <= 0.0 || alpha >= 1.0 {
1358        return None;
1359    }
1360
1361    // Fit original model
1362    let original = functional_logistic(data, y, scalar_covariates, ncomp, max_iter, tol)?;
1363    let center = original.beta_t.clone();
1364
1365    // Bootstrap replicates
1366    let boot_betas: Vec<Vec<f64>> = iter_maybe_parallel!(0..n_boot)
1367        .filter_map(|b| {
1368            let mut rng = StdRng::seed_from_u64(seed.wrapping_add(b as u64));
1369            let indices: Vec<usize> = (0..n).map(|_| rng.gen_range(0..n)).collect();
1370
1371            let boot_data = subsample_rows(data, &indices);
1372            let boot_y: Vec<f64> = indices.iter().map(|&i| y[i]).collect();
1373            let boot_sc = scalar_covariates.map(|sc| subsample_rows(sc, &indices));
1374
1375            functional_logistic(&boot_data, &boot_y, boot_sc.as_ref(), ncomp, max_iter, tol)
1376                .map(|fit| fit.beta_t)
1377        })
1378        .collect();
1379
1380    let n_boot_success = boot_betas.len();
1381    if n_boot_success < 3 {
1382        return None;
1383    }
1384
1385    let lo_q = alpha / 2.0;
1386    let hi_q = 1.0 - alpha / 2.0;
1387    let mut lower = vec![0.0; m];
1388    let mut upper = vec![0.0; m];
1389    let mut boot_se = vec![0.0; m];
1390
1391    for j in 0..m {
1392        let mut vals: Vec<f64> = boot_betas.iter().map(|b| b[j]).collect();
1393        vals.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
1394        lower[j] = quantile(&vals, lo_q);
1395        upper[j] = quantile(&vals, hi_q);
1396
1397        let mean_j: f64 = vals.iter().sum::<f64>() / vals.len() as f64;
1398        let var_j: f64 =
1399            vals.iter().map(|&v| (v - mean_j).powi(2)).sum::<f64>() / (vals.len() - 1) as f64;
1400        boot_se[j] = var_j.sqrt().max(1e-15);
1401    }
1402
1403    let mut t_stats: Vec<f64> = boot_betas
1404        .iter()
1405        .map(|b| {
1406            (0..m)
1407                .map(|j| ((b[j] - center[j]) / boot_se[j]).abs())
1408                .fold(0.0_f64, f64::max)
1409        })
1410        .collect();
1411    t_stats.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
1412
1413    let c_alpha = quantile(&t_stats, 1.0 - alpha);
1414    let sim_lower: Vec<f64> = (0..m).map(|j| center[j] - c_alpha * boot_se[j]).collect();
1415    let sim_upper: Vec<f64> = (0..m).map(|j| center[j] + c_alpha * boot_se[j]).collect();
1416
1417    Some(BootstrapCiResult {
1418        lower,
1419        upper,
1420        center,
1421        sim_lower,
1422        sim_upper,
1423        n_boot_success,
1424    })
1425}
1426
1427// ---------------------------------------------------------------------------
1428// Basis regression CV (R's fregre.basis.cv)
1429// ---------------------------------------------------------------------------
1430
1431/// Result of lambda selection for basis regression via cross-validation.
1432pub struct FregreBasisCvResult {
1433    /// Optimal smoothing parameter lambda.
1434    pub optimal_lambda: f64,
1435    /// Mean CV error for each lambda.
1436    pub cv_errors: Vec<f64>,
1437    /// SE of CV error across folds for each lambda.
1438    pub cv_se: Vec<f64>,
1439    /// Lambda values tested.
1440    pub lambda_values: Vec<f64>,
1441    /// Minimum mean CV error.
1442    pub min_cv_error: f64,
1443}
1444
1445/// K-fold CV for selecting the regularization parameter lambda
1446/// in basis-regression (R's `fregre.basis.cv`).
1447///
1448/// Projects functional predictors onto a B-spline or Fourier basis,
1449/// fits a ridge regression on the resulting coefficients for each lambda,
1450/// and selects the lambda minimizing mean CV error.
1451///
1452/// # Arguments
1453/// * `data` — Functional predictor matrix (n × m)
1454/// * `y` — Scalar response vector (length n)
1455/// * `argvals` — Evaluation grid (length m)
1456/// * `n_folds` — Number of CV folds
1457/// * `lambda_range` — Lambda values to test. If `None`, uses 20 log-spaced
1458///   values from 1e-4 to 1e4.
1459/// * `nbasis` — Number of basis functions
1460/// * `basis_type` — Basis type (B-spline or Fourier)
1461pub fn fregre_basis_cv(
1462    data: &FdMatrix,
1463    y: &[f64],
1464    argvals: &[f64],
1465    n_folds: usize,
1466    lambda_range: Option<&[f64]>,
1467    nbasis: usize,
1468    basis_type: &crate::smooth_basis::BasisType,
1469) -> Option<FregreBasisCvResult> {
1470    use crate::smooth_basis::{smooth_basis, BasisType, FdPar};
1471
1472    let (n, m) = data.shape();
1473    if n < n_folds || m == 0 || y.len() != n || argvals.len() != m || nbasis < 2 {
1474        return None;
1475    }
1476
1477    // Default lambda grid
1478    let default_lambdas: Vec<f64> = if lambda_range.is_none() {
1479        (0..20)
1480            .map(|i| {
1481                let log_lam = -4.0 + 8.0 * i as f64 / 19.0;
1482                10.0_f64.powf(log_lam)
1483            })
1484            .collect()
1485    } else {
1486        Vec::new()
1487    };
1488    let lambdas = lambda_range.unwrap_or(&default_lambdas);
1489
1490    if lambdas.is_empty() {
1491        return None;
1492    }
1493
1494    // Compute penalty matrix once
1495    let penalty = match basis_type {
1496        BasisType::Bspline { order } => {
1497            crate::smooth_basis::bspline_penalty_matrix(argvals, nbasis, *order, 2)
1498        }
1499        BasisType::Fourier { period } => {
1500            crate::smooth_basis::fourier_penalty_matrix(nbasis, *period, 2)
1501        }
1502    };
1503
1504    // Create folds
1505    let folds = crate::cv::create_folds(n, n_folds, 42);
1506
1507    let mut cv_errors = vec![0.0; lambdas.len()];
1508    let mut cv_fold_errors: Vec<Vec<f64>> = vec![Vec::with_capacity(n_folds); lambdas.len()];
1509
1510    for fold in 0..n_folds {
1511        let (train_idx, test_idx) = crate::cv::fold_indices(&folds, fold);
1512        if train_idx.is_empty() || test_idx.is_empty() {
1513            continue;
1514        }
1515
1516        let train_data = crate::cv::subset_rows(data, &train_idx);
1517        let train_y = crate::cv::subset_vec(y, &train_idx);
1518        let test_data = crate::cv::subset_rows(data, &test_idx);
1519        let test_y = crate::cv::subset_vec(y, &test_idx);
1520
1521        for (li, &lam) in lambdas.iter().enumerate() {
1522            let fdpar = FdPar {
1523                basis_type: basis_type.clone(),
1524                nbasis,
1525                lambda: lam,
1526                lfd_order: 2,
1527                penalty_matrix: penalty.clone(),
1528            };
1529
1530            // Project training curves to basis coefficients
1531            let train_result = match smooth_basis(&train_data, argvals, &fdpar) {
1532                Some(r) => r,
1533                None => {
1534                    cv_fold_errors[li].push(f64::INFINITY);
1535                    continue;
1536                }
1537            };
1538
1539            // Project test curves to basis coefficients
1540            let test_result = match smooth_basis(&test_data, argvals, &fdpar) {
1541                Some(r) => r,
1542                None => {
1543                    cv_fold_errors[li].push(f64::INFINITY);
1544                    continue;
1545                }
1546            };
1547
1548            // Fit linear model on training coefficients -> y
1549            let train_coefs = &train_result.coefficients;
1550            let test_coefs = &test_result.coefficients;
1551            let n_train = train_idx.len();
1552            let n_test = test_idx.len();
1553            let k = train_coefs.ncols();
1554
1555            // Penalized OLS: (X'X + lam*P) \ X'y
1556            let mut xtx = vec![0.0; k * k];
1557            let mut xty_vec = vec![0.0; k];
1558            for i in 0..n_train {
1559                for j in 0..k {
1560                    xty_vec[j] += train_coefs[(i, j)] * train_y[i];
1561                    for l in 0..k {
1562                        xtx[j * k + l] += train_coefs[(i, j)] * train_coefs[(i, l)];
1563                    }
1564                }
1565            }
1566            // Add roughness penalty: lam * P (not ridge lam * I)
1567            for j in 0..k {
1568                for l in 0..k {
1569                    xtx[j * k + l] += lam * penalty[j * k + l];
1570                }
1571            }
1572
1573            // Solve via Gaussian elimination
1574            let beta = {
1575                let mut a = xtx;
1576                let mut b = xty_vec;
1577                crate::smoothing::solve_gaussian_pub(&mut a, &mut b, k)
1578            };
1579
1580            // Predict on test set
1581            let mut fold_mse = 0.0;
1582            for i in 0..n_test {
1583                let mut yhat = 0.0;
1584                for j in 0..k {
1585                    yhat += test_coefs[(i, j)] * beta[j];
1586                }
1587                fold_mse += (test_y[i] - yhat).powi(2);
1588            }
1589            fold_mse /= n_test as f64;
1590            cv_fold_errors[li].push(fold_mse);
1591        }
1592    }
1593
1594    // Compute mean and SE across folds
1595    let mut cv_se = vec![0.0; lambdas.len()];
1596    for li in 0..lambdas.len() {
1597        let errs = &cv_fold_errors[li];
1598        if errs.is_empty() {
1599            cv_errors[li] = f64::INFINITY;
1600            continue;
1601        }
1602        let mean = errs.iter().sum::<f64>() / errs.len() as f64;
1603        cv_errors[li] = mean;
1604        if errs.len() > 1 {
1605            let var =
1606                errs.iter().map(|&e| (e - mean).powi(2)).sum::<f64>() / (errs.len() - 1) as f64;
1607            cv_se[li] = (var / errs.len() as f64).sqrt();
1608        }
1609    }
1610
1611    let (best_idx, &min_cv_error) = cv_errors
1612        .iter()
1613        .enumerate()
1614        .min_by(|(_, a), (_, b)| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal))?;
1615
1616    Some(FregreBasisCvResult {
1617        optimal_lambda: lambdas[best_idx],
1618        cv_errors,
1619        cv_se,
1620        lambda_values: lambdas.to_vec(),
1621        min_cv_error,
1622    })
1623}
1624
1625// ---------------------------------------------------------------------------
1626// Nonparametric regression bandwidth CV (R's fregre.np.cv)
1627// ---------------------------------------------------------------------------
1628
1629/// Result of bandwidth selection for nonparametric regression via CV.
1630pub struct FregreNpCvResult {
1631    /// Optimal bandwidth.
1632    pub optimal_h: f64,
1633    /// Mean CV error for each bandwidth.
1634    pub cv_errors: Vec<f64>,
1635    /// SE of CV error across folds for each bandwidth.
1636    pub cv_se: Vec<f64>,
1637    /// Bandwidth values tested.
1638    pub h_values: Vec<f64>,
1639    /// Minimum mean CV error.
1640    pub min_cv_error: f64,
1641}
1642
1643/// K-fold CV for selecting the bandwidth in functional nonparametric
1644/// regression (R's `fregre.np.cv`).
1645///
1646/// Computes a full n×n L2 distance matrix once, then for each candidate
1647/// bandwidth and each fold, does Nadaraya-Watson prediction.
1648///
1649/// # Arguments
1650/// * `data` — Functional predictor matrix (n × m)
1651/// * `y` — Scalar response vector (length n)
1652/// * `argvals` — Evaluation grid (length m)
1653/// * `n_folds` — Number of CV folds
1654/// * `h_range` — Bandwidth values to test. If `None`, uses 20 quantiles
1655///   (5th to 95th percentile) of the pairwise L2 distances.
1656/// * `scalar_covariates` — Optional scalar covariates (n × p)
1657pub fn fregre_np_cv(
1658    data: &FdMatrix,
1659    y: &[f64],
1660    argvals: &[f64],
1661    n_folds: usize,
1662    h_range: Option<&[f64]>,
1663    scalar_covariates: Option<&FdMatrix>,
1664) -> Option<FregreNpCvResult> {
1665    let (n, m) = data.shape();
1666    if n < n_folds || m == 0 || y.len() != n || argvals.len() != m || n < 3 {
1667        return None;
1668    }
1669
1670    // Compute full distance matrix
1671    let func_dists = compute_pairwise_distances(data, argvals);
1672    let has_scalar = scalar_covariates.is_some();
1673    let scalar_dists = scalar_covariates
1674        .map(compute_scalar_distances)
1675        .unwrap_or_default();
1676
1677    // Default h grid: 20 quantiles of pairwise distances
1678    let default_h: Vec<f64> = if h_range.is_none() {
1679        let mut nonzero: Vec<f64> = Vec::new();
1680        for i in 0..n {
1681            for j in (i + 1)..n {
1682                let d = func_dists[i * n + j];
1683                if d > 0.0 {
1684                    nonzero.push(d);
1685                }
1686            }
1687        }
1688        nonzero.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
1689        if nonzero.is_empty() {
1690            return None;
1691        }
1692        (1..=20)
1693            .map(|qi| {
1694                let q = 0.05 + 0.90 * (qi - 1) as f64 / 19.0;
1695                let idx = ((nonzero.len() as f64 * q) as usize).min(nonzero.len() - 1);
1696                nonzero[idx].max(1e-10)
1697            })
1698            .collect()
1699    } else {
1700        Vec::new()
1701    };
1702    let h_values = h_range.unwrap_or(&default_h);
1703
1704    if h_values.is_empty() {
1705        return None;
1706    }
1707
1708    let folds = crate::cv::create_folds(n, n_folds, 42);
1709
1710    let mut cv_errors = vec![0.0; h_values.len()];
1711    let mut cv_fold_errors: Vec<Vec<f64>> = vec![Vec::with_capacity(n_folds); h_values.len()];
1712
1713    // Select scalar bandwidth once (from full data)
1714    let h_scalar = if has_scalar {
1715        select_bandwidth_loo(&scalar_dists, y, n, Some(&func_dists))
1716    } else {
1717        1.0
1718    };
1719
1720    for fold in 0..n_folds {
1721        let (train_idx, test_idx) = crate::cv::fold_indices(&folds, fold);
1722        if train_idx.is_empty() || test_idx.is_empty() {
1723            continue;
1724        }
1725
1726        let train_y: Vec<f64> = train_idx.iter().map(|&i| y[i]).collect();
1727
1728        for (hi, &h) in h_values.iter().enumerate() {
1729            let mut fold_mse = 0.0;
1730            for &ti in &test_idx {
1731                // NW prediction using only training data
1732                let mut num = 0.0;
1733                let mut den = 0.0;
1734                for (local_j, &tj) in train_idx.iter().enumerate() {
1735                    let kf = gaussian_kernel(func_dists[ti * n + tj], h);
1736                    let ks = if has_scalar {
1737                        gaussian_kernel(scalar_dists[ti * n + tj], h_scalar)
1738                    } else {
1739                        1.0
1740                    };
1741                    let w = kf * ks;
1742                    num += w * train_y[local_j];
1743                    den += w;
1744                }
1745                let y_hat = if den > 1e-15 { num / den } else { y[ti] };
1746                fold_mse += (y[ti] - y_hat).powi(2);
1747            }
1748            fold_mse /= test_idx.len() as f64;
1749            cv_fold_errors[hi].push(fold_mse);
1750        }
1751    }
1752
1753    // Compute mean and SE
1754    let mut cv_se = vec![0.0; h_values.len()];
1755    for hi in 0..h_values.len() {
1756        let errs = &cv_fold_errors[hi];
1757        if errs.is_empty() {
1758            cv_errors[hi] = f64::INFINITY;
1759            continue;
1760        }
1761        let mean = errs.iter().sum::<f64>() / errs.len() as f64;
1762        cv_errors[hi] = mean;
1763        if errs.len() > 1 {
1764            let var =
1765                errs.iter().map(|&e| (e - mean).powi(2)).sum::<f64>() / (errs.len() - 1) as f64;
1766            cv_se[hi] = (var / errs.len() as f64).sqrt();
1767        }
1768    }
1769
1770    let (best_idx, &min_cv_error) = cv_errors
1771        .iter()
1772        .enumerate()
1773        .min_by(|(_, a), (_, b)| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal))?;
1774
1775    Some(FregreNpCvResult {
1776        optimal_h: h_values[best_idx],
1777        cv_errors,
1778        cv_se,
1779        h_values: h_values.to_vec(),
1780        min_cv_error,
1781    })
1782}
1783
1784/// Predict probabilities P(Y=1) for new data using a fitted functional logistic model.
1785///
1786/// Projects new curves through the stored FPCA, computes linear predictor,
1787/// and applies sigmoid.
1788///
1789/// # Arguments
1790/// * `fit` — A fitted [`FunctionalLogisticResult`]
1791/// * `new_data` — New functional predictor matrix (n_new × m)
1792/// * `new_scalar` — Optional new scalar covariates (n_new × p)
1793pub fn predict_functional_logistic(
1794    fit: &FunctionalLogisticResult,
1795    new_data: &FdMatrix,
1796    new_scalar: Option<&FdMatrix>,
1797) -> Vec<f64> {
1798    let (n_new, m) = new_data.shape();
1799    let ncomp = fit.ncomp;
1800    let p_scalar = fit.gamma.len();
1801
1802    (0..n_new)
1803        .map(|i| {
1804            let mut eta = fit.coefficients[0]; // intercept
1805            for k in 0..ncomp {
1806                let mut s = 0.0;
1807                for j in 0..m {
1808                    s += (new_data[(i, j)] - fit.fpca.mean[j]) * fit.fpca.rotation[(j, k)];
1809                }
1810                eta += fit.coefficients[1 + k] * s;
1811            }
1812            if let Some(sc) = new_scalar {
1813                for j in 0..p_scalar {
1814                    eta += fit.gamma[j] * sc[(i, j)];
1815                }
1816            }
1817            sigmoid(eta)
1818        })
1819        .collect()
1820}
1821
1822// ---------------------------------------------------------------------------
1823// Predict methods on result structs
1824// ---------------------------------------------------------------------------
1825
1826impl FregreLmResult {
1827    /// Predict new responses. Delegates to [`predict_fregre_lm`].
1828    pub fn predict(&self, new_data: &FdMatrix, new_scalar: Option<&FdMatrix>) -> Vec<f64> {
1829        predict_fregre_lm(self, new_data, new_scalar)
1830    }
1831}
1832
1833impl FunctionalLogisticResult {
1834    /// Predict P(Y=1) for new data. Delegates to [`predict_functional_logistic`].
1835    pub fn predict(&self, new_data: &FdMatrix, new_scalar: Option<&FdMatrix>) -> Vec<f64> {
1836        predict_functional_logistic(self, new_data, new_scalar)
1837    }
1838}
1839
1840// ---------------------------------------------------------------------------
1841// Tests
1842// ---------------------------------------------------------------------------
1843
1844#[cfg(test)]
1845mod tests {
1846    use super::*;
1847    use std::f64::consts::PI;
1848
1849    fn generate_test_data(n: usize, m: usize, seed: u64) -> (FdMatrix, Vec<f64>, Vec<f64>) {
1850        let t: Vec<f64> = (0..m).map(|j| j as f64 / (m - 1) as f64).collect();
1851        let mut data = FdMatrix::zeros(n, m);
1852        let mut y = vec![0.0; n];
1853
1854        for i in 0..n {
1855            let phase =
1856                (seed.wrapping_mul(17).wrapping_add(i as u64 * 31) % 1000) as f64 / 1000.0 * PI;
1857            let amplitude =
1858                ((seed.wrapping_mul(13).wrapping_add(i as u64 * 7) % 100) as f64 / 100.0) - 0.5;
1859            for j in 0..m {
1860                data[(i, j)] =
1861                    (2.0 * PI * t[j] + phase).sin() + amplitude * (4.0 * PI * t[j]).cos();
1862            }
1863            y[i] = 2.0 * phase + 3.0 * amplitude + 0.05 * (seed.wrapping_add(i as u64) % 10) as f64;
1864        }
1865        (data, y, t)
1866    }
1867
1868    // ----- fregre_lm tests -----
1869
1870    #[test]
1871    fn test_fregre_lm_basic() {
1872        let (data, y, _t) = generate_test_data(30, 50, 42);
1873        let result = fregre_lm(&data, &y, None, 3);
1874        assert!(result.is_some());
1875        let fit = result.unwrap();
1876        assert_eq!(fit.fitted_values.len(), 30);
1877        assert_eq!(fit.residuals.len(), 30);
1878        assert_eq!(fit.beta_t.len(), 50);
1879        assert_eq!(fit.ncomp, 3);
1880        assert!(fit.r_squared >= 0.0 && fit.r_squared <= 1.0 + 1e-10);
1881    }
1882
1883    #[test]
1884    fn test_fregre_lm_with_scalar_covariates() {
1885        let (data, y, _t) = generate_test_data(30, 50, 42);
1886        let mut sc = FdMatrix::zeros(30, 2);
1887        for i in 0..30 {
1888            sc[(i, 0)] = i as f64 / 30.0;
1889            sc[(i, 1)] = (i as f64 * 0.7).sin();
1890        }
1891        let result = fregre_lm(&data, &y, Some(&sc), 3);
1892        assert!(result.is_some());
1893        let fit = result.unwrap();
1894        assert_eq!(fit.gamma.len(), 2);
1895    }
1896
1897    #[test]
1898    fn test_fregre_lm_residuals_sum_near_zero() {
1899        let (data, y, _t) = generate_test_data(30, 50, 42);
1900        let fit = fregre_lm(&data, &y, None, 3).unwrap();
1901        let resid_sum: f64 = fit.residuals.iter().sum::<f64>();
1902        assert!(
1903            resid_sum.abs() < 1e-8,
1904            "Residuals should sum to ~0 with intercept, got {}",
1905            resid_sum
1906        );
1907    }
1908
1909    #[test]
1910    fn test_fregre_lm_fitted_plus_residuals_equals_y() {
1911        let (data, y, _t) = generate_test_data(30, 50, 42);
1912        let fit = fregre_lm(&data, &y, None, 3).unwrap();
1913        for i in 0..30 {
1914            let reconstructed = fit.fitted_values[i] + fit.residuals[i];
1915            assert!(
1916                (reconstructed - y[i]).abs() < 1e-10,
1917                "ŷ + r should equal y at index {}",
1918                i
1919            );
1920        }
1921    }
1922
1923    #[test]
1924    fn test_fregre_lm_more_components_higher_r2() {
1925        let (data, y, _t) = generate_test_data(50, 50, 42);
1926        let fit1 = fregre_lm(&data, &y, None, 1).unwrap();
1927        let fit3 = fregre_lm(&data, &y, None, 3).unwrap();
1928        assert!(
1929            fit3.r_squared >= fit1.r_squared - 1e-10,
1930            "More components should give >= R²: {} vs {}",
1931            fit3.r_squared,
1932            fit1.r_squared
1933        );
1934    }
1935
1936    #[test]
1937    fn test_fregre_lm_invalid_input() {
1938        let data = FdMatrix::zeros(2, 50);
1939        let y = vec![1.0, 2.0];
1940        assert!(fregre_lm(&data, &y, None, 1).is_none());
1941
1942        let data = FdMatrix::zeros(10, 50);
1943        let y = vec![1.0; 5];
1944        assert!(fregre_lm(&data, &y, None, 2).is_none());
1945    }
1946
1947    #[test]
1948    fn test_fregre_lm_std_errors_positive() {
1949        let (data, y, _t) = generate_test_data(30, 50, 42);
1950        let fit = fregre_lm(&data, &y, None, 3).unwrap();
1951        for (i, &se) in fit.std_errors.iter().enumerate() {
1952            assert!(
1953                se > 0.0 && se.is_finite(),
1954                "Std error {} should be positive finite, got {}",
1955                i,
1956                se
1957            );
1958        }
1959    }
1960
1961    // ----- predict_fregre_lm tests -----
1962
1963    #[test]
1964    fn test_predict_fregre_lm_on_training_data() {
1965        let (data, y, _t) = generate_test_data(30, 50, 42);
1966        let fit = fregre_lm(&data, &y, None, 3).unwrap();
1967        let preds = predict_fregre_lm(&fit, &data, None);
1968        for i in 0..30 {
1969            assert!(
1970                (preds[i] - fit.fitted_values[i]).abs() < 1e-6,
1971                "Prediction on training data should match fitted values"
1972            );
1973        }
1974    }
1975
1976    // ----- fregre_cv tests -----
1977
1978    #[test]
1979    fn test_fregre_cv_returns_result() {
1980        let (data, y, _t) = generate_test_data(30, 50, 42);
1981        let result = fregre_cv(&data, &y, None, 1, 8, 5);
1982        assert!(result.is_some());
1983        let cv = result.unwrap();
1984        assert!(cv.optimal_k >= 1 && cv.optimal_k <= 8);
1985        assert!(cv.min_cv_error >= 0.0);
1986    }
1987
1988    #[test]
1989    fn test_fregre_cv_with_scalar_covariates() {
1990        let (data, y, _t) = generate_test_data(30, 50, 42);
1991        let mut sc = FdMatrix::zeros(30, 1);
1992        for i in 0..30 {
1993            sc[(i, 0)] = i as f64;
1994        }
1995        let result = fregre_cv(&data, &y, Some(&sc), 1, 5, 3);
1996        assert!(result.is_some());
1997    }
1998
1999    // ----- fregre_np_mixed tests -----
2000
2001    #[test]
2002    fn test_fregre_np_mixed_basic() {
2003        let (data, y, t) = generate_test_data(30, 50, 42);
2004        let result = fregre_np_mixed(&data, &y, &t, None, 0.0, 0.0);
2005        assert!(result.is_some());
2006        let fit = result.unwrap();
2007        assert_eq!(fit.fitted_values.len(), 30);
2008        assert!(fit.h_func > 0.0);
2009        assert!(fit.cv_error >= 0.0);
2010    }
2011
2012    #[test]
2013    fn test_fregre_np_mixed_with_scalars() {
2014        let (data, y, t) = generate_test_data(30, 50, 42);
2015        let mut sc = FdMatrix::zeros(30, 1);
2016        for i in 0..30 {
2017            sc[(i, 0)] = i as f64 / 30.0;
2018        }
2019        let result = fregre_np_mixed(&data, &y, &t, Some(&sc), 0.0, 0.0);
2020        assert!(result.is_some());
2021        let fit = result.unwrap();
2022        assert!(fit.h_scalar > 0.0);
2023    }
2024
2025    #[test]
2026    fn test_fregre_np_mixed_manual_bandwidth() {
2027        let (data, y, t) = generate_test_data(30, 50, 42);
2028        let result = fregre_np_mixed(&data, &y, &t, None, 0.5, 0.0);
2029        assert!(result.is_some());
2030        let fit = result.unwrap();
2031        assert!(
2032            (fit.h_func - 0.5).abs() < 1e-10,
2033            "Should use provided bandwidth"
2034        );
2035    }
2036
2037    // ----- functional_logistic tests -----
2038
2039    #[test]
2040    fn test_functional_logistic_basic() {
2041        let (data, y_cont, _t) = generate_test_data(30, 50, 42);
2042        let y_median = {
2043            let mut sorted = y_cont.clone();
2044            sorted.sort_by(|a, b| a.partial_cmp(b).unwrap());
2045            sorted[sorted.len() / 2]
2046        };
2047        let y_bin: Vec<f64> = y_cont
2048            .iter()
2049            .map(|&v| if v >= y_median { 1.0 } else { 0.0 })
2050            .collect();
2051
2052        let result = functional_logistic(&data, &y_bin, None, 3, 25, 1e-6);
2053        assert!(result.is_some());
2054        let fit = result.unwrap();
2055        assert_eq!(fit.probabilities.len(), 30);
2056        assert_eq!(fit.predicted_classes.len(), 30);
2057        assert!(fit.accuracy >= 0.0 && fit.accuracy <= 1.0);
2058        for &p in &fit.probabilities {
2059            assert!((0.0..=1.0).contains(&p), "Probability out of range: {}", p);
2060        }
2061    }
2062
2063    #[test]
2064    fn test_functional_logistic_with_scalar_covariates() {
2065        let (data, y_cont, _t) = generate_test_data(30, 50, 42);
2066        let y_median = {
2067            let mut sorted = y_cont.clone();
2068            sorted.sort_by(|a, b| a.partial_cmp(b).unwrap());
2069            sorted[sorted.len() / 2]
2070        };
2071        let y_bin: Vec<f64> = y_cont
2072            .iter()
2073            .map(|&v| if v >= y_median { 1.0 } else { 0.0 })
2074            .collect();
2075
2076        let mut sc = FdMatrix::zeros(30, 1);
2077        for i in 0..30 {
2078            sc[(i, 0)] = i as f64 / 30.0;
2079        }
2080        let result = functional_logistic(&data, &y_bin, Some(&sc), 3, 25, 1e-6);
2081        assert!(result.is_some());
2082        let fit = result.unwrap();
2083        assert_eq!(fit.gamma.len(), 1);
2084    }
2085
2086    #[test]
2087    fn test_functional_logistic_invalid_response() {
2088        let (data, _, _) = generate_test_data(30, 50, 42);
2089        let y: Vec<f64> = (0..30).map(|i| i as f64).collect();
2090        assert!(functional_logistic(&data, &y, None, 3, 25, 1e-6).is_none());
2091    }
2092
2093    #[test]
2094    fn test_functional_logistic_convergence() {
2095        let (data, y_cont, _t) = generate_test_data(40, 50, 42);
2096        let y_median = {
2097            let mut sorted = y_cont.clone();
2098            sorted.sort_by(|a, b| a.partial_cmp(b).unwrap());
2099            sorted[sorted.len() / 2]
2100        };
2101        let y_bin: Vec<f64> = y_cont
2102            .iter()
2103            .map(|&v| if v >= y_median { 1.0 } else { 0.0 })
2104            .collect();
2105
2106        let fit = functional_logistic(&data, &y_bin, None, 3, 100, 1e-8).unwrap();
2107        assert!(fit.iterations <= 100, "Should converge within max_iter");
2108    }
2109
2110    // ----- Edge cases -----
2111
2112    #[test]
2113    fn test_fregre_lm_single_component() {
2114        let (data, y, _t) = generate_test_data(20, 50, 42);
2115        let result = fregre_lm(&data, &y, None, 1);
2116        assert!(result.is_some());
2117        let fit = result.unwrap();
2118        assert_eq!(fit.ncomp, 1);
2119    }
2120
2121    #[test]
2122    fn test_fregre_lm_auto_k_selection() {
2123        let (data, y, _t) = generate_test_data(30, 50, 42);
2124        let result = fregre_lm(&data, &y, None, 0);
2125        assert!(result.is_some());
2126        let fit = result.unwrap();
2127        assert!(fit.ncomp >= 1);
2128    }
2129
2130    #[test]
2131    fn test_predict_fregre_np_basic() {
2132        let (data, y, t) = generate_test_data(30, 50, 42);
2133        let preds = predict_fregre_np(&data, &y, None, &data, None, &t, 0.5, 0.5);
2134        assert_eq!(preds.len(), 30);
2135        for &p in &preds {
2136            assert!(p.is_finite(), "Prediction should be finite");
2137        }
2138    }
2139
2140    #[test]
2141    fn test_sigmoid_properties() {
2142        assert!((sigmoid(0.0) - 0.5).abs() < 1e-10);
2143        assert!(sigmoid(10.0) > 0.999);
2144        assert!(sigmoid(-10.0) < 0.001);
2145        assert!((sigmoid(3.0) + sigmoid(-3.0) - 1.0).abs() < 1e-10);
2146    }
2147
2148    // ----- beta_se tests -----
2149
2150    #[test]
2151    fn test_fregre_lm_beta_se() {
2152        let (data, y, _t) = generate_test_data(30, 50, 42);
2153        let fit = fregre_lm(&data, &y, None, 3).unwrap();
2154        assert_eq!(fit.beta_se.len(), 50, "beta_se should have length m");
2155        for (j, &se) in fit.beta_se.iter().enumerate() {
2156            assert!(
2157                se > 0.0 && se.is_finite(),
2158                "beta_se[{}] should be positive finite, got {}",
2159                j,
2160                se
2161            );
2162        }
2163    }
2164
2165    #[test]
2166    fn test_fregre_lm_beta_se_coverage() {
2167        // Use generate_test_data which is known to produce valid FPCA, then check SE properties
2168        let (data, y, _t) = generate_test_data(50, 50, 99);
2169        let fit = fregre_lm(&data, &y, None, 3).unwrap();
2170        assert_eq!(fit.beta_se.len(), 50);
2171        // With valid data, beta_se should be positive everywhere
2172        for (j, &se) in fit.beta_se.iter().enumerate() {
2173            assert!(
2174                se > 0.0 && se.is_finite(),
2175                "beta_se[{}] should be positive finite, got {}",
2176                j,
2177                se
2178            );
2179        }
2180        // The CI band [beta_t ± 2·SE] should have non-zero width everywhere
2181        for j in 0..50 {
2182            let width = 4.0 * fit.beta_se[j];
2183            assert!(width > 0.0, "CI width should be positive at j={}", j);
2184        }
2185    }
2186
2187    #[test]
2188    fn test_functional_logistic_beta_se() {
2189        let (data, y_cont, _t) = generate_test_data(30, 50, 42);
2190        let y_median = {
2191            let mut sorted = y_cont.clone();
2192            sorted.sort_by(|a, b| a.partial_cmp(b).unwrap());
2193            sorted[sorted.len() / 2]
2194        };
2195        let y_bin: Vec<f64> = y_cont
2196            .iter()
2197            .map(|&v| if v >= y_median { 1.0 } else { 0.0 })
2198            .collect();
2199
2200        let fit = functional_logistic(&data, &y_bin, None, 3, 25, 1e-6).unwrap();
2201        assert_eq!(fit.beta_se.len(), 50, "beta_se should have length m");
2202        assert_eq!(
2203            fit.std_errors.len(),
2204            1 + 3,
2205            "std_errors should have length 1 + ncomp"
2206        );
2207        for (j, &se) in fit.beta_se.iter().enumerate() {
2208            assert!(
2209                se > 0.0 && se.is_finite(),
2210                "beta_se[{}] should be positive finite, got {}",
2211                j,
2212                se
2213            );
2214        }
2215        for (j, &se) in fit.std_errors.iter().enumerate() {
2216            assert!(
2217                se > 0.0 && se.is_finite(),
2218                "std_errors[{}] should be positive finite, got {}",
2219                j,
2220                se
2221            );
2222        }
2223    }
2224
2225    #[test]
2226    fn test_beta_se_zero_for_constant() {
2227        // When all curves are nearly identical, β(t) ≈ 0 but SE should still be finite/positive
2228        let n = 30;
2229        let m = 20;
2230        let mut data = FdMatrix::zeros(n, m);
2231        let mut y = vec![0.0; n];
2232        for i in 0..n {
2233            for j in 0..m {
2234                // Nearly identical curves with tiny variation
2235                data[(i, j)] = 1.0 + 0.001 * (i as f64 / n as f64);
2236            }
2237            y[i] = i as f64 / n as f64;
2238        }
2239        let fit = fregre_lm(&data, &y, None, 1).unwrap();
2240        assert_eq!(fit.beta_se.len(), m);
2241        for (j, &se) in fit.beta_se.iter().enumerate() {
2242            assert!(
2243                se.is_finite() && se >= 0.0,
2244                "beta_se[{}] should be finite non-negative, got {}",
2245                j,
2246                se
2247            );
2248        }
2249    }
2250
2251    // ----- Bootstrap CI tests -----
2252
2253    #[test]
2254    fn test_bootstrap_ci_fregre_lm_shape() {
2255        let (data, y, _t) = generate_test_data(30, 20, 42);
2256        let result = bootstrap_ci_fregre_lm(&data, &y, None, 2, 50, 0.05, 123);
2257        assert!(result.is_some(), "bootstrap_ci_fregre_lm should succeed");
2258        let ci = result.unwrap();
2259        assert_eq!(ci.lower.len(), 20);
2260        assert_eq!(ci.upper.len(), 20);
2261        assert_eq!(ci.center.len(), 20);
2262        assert_eq!(ci.sim_lower.len(), 20);
2263        assert_eq!(ci.sim_upper.len(), 20);
2264        assert!(ci.n_boot_success > 0);
2265    }
2266
2267    #[test]
2268    fn test_bootstrap_ci_fregre_lm_ordering() {
2269        let (data, y, _t) = generate_test_data(30, 20, 42);
2270        let ci = bootstrap_ci_fregre_lm(&data, &y, None, 2, 100, 0.05, 42).unwrap();
2271        for j in 0..20 {
2272            // Pointwise: lower ≤ center ≤ upper
2273            assert!(
2274                ci.lower[j] <= ci.center[j] + 1e-10,
2275                "lower <= center at j={}: {} > {}",
2276                j,
2277                ci.lower[j],
2278                ci.center[j]
2279            );
2280            assert!(
2281                ci.center[j] <= ci.upper[j] + 1e-10,
2282                "center <= upper at j={}: {} > {}",
2283                j,
2284                ci.center[j],
2285                ci.upper[j]
2286            );
2287            // Simultaneous: sim_lower ≤ center ≤ sim_upper
2288            assert!(
2289                ci.sim_lower[j] <= ci.center[j] + 1e-10,
2290                "sim_lower <= center at j={}: {} > {}",
2291                j,
2292                ci.sim_lower[j],
2293                ci.center[j]
2294            );
2295            assert!(
2296                ci.center[j] <= ci.sim_upper[j] + 1e-10,
2297                "center <= sim_upper at j={}: {} > {}",
2298                j,
2299                ci.center[j],
2300                ci.sim_upper[j]
2301            );
2302        }
2303        // Simultaneous band should be wider on average
2304        let pw_width: f64 = (0..20).map(|j| ci.upper[j] - ci.lower[j]).sum::<f64>() / 20.0;
2305        let sim_width: f64 = (0..20)
2306            .map(|j| ci.sim_upper[j] - ci.sim_lower[j])
2307            .sum::<f64>()
2308            / 20.0;
2309        assert!(
2310            sim_width >= pw_width - 1e-10,
2311            "Simultaneous band should be wider on average: sim={}, pw={}",
2312            sim_width,
2313            pw_width
2314        );
2315    }
2316
2317    #[test]
2318    fn test_bootstrap_ci_fregre_lm_center_matches_fit() {
2319        let (data, y, _t) = generate_test_data(30, 20, 42);
2320        let fit = fregre_lm(&data, &y, None, 2).unwrap();
2321        let ci = bootstrap_ci_fregre_lm(&data, &y, None, 2, 50, 0.05, 42).unwrap();
2322        for j in 0..20 {
2323            assert!(
2324                (ci.center[j] - fit.beta_t[j]).abs() < 1e-12,
2325                "center should match original beta_t at j={}",
2326                j
2327            );
2328        }
2329    }
2330
2331    #[test]
2332    fn test_bootstrap_ci_functional_logistic_shape() {
2333        let (data, y_cont, _t) = generate_test_data(40, 20, 42);
2334        let y_median = {
2335            let mut sorted = y_cont.clone();
2336            sorted.sort_by(|a, b| a.partial_cmp(b).unwrap());
2337            sorted[sorted.len() / 2]
2338        };
2339        let y_bin: Vec<f64> = y_cont
2340            .iter()
2341            .map(|&v| if v >= y_median { 1.0 } else { 0.0 })
2342            .collect();
2343
2344        let result =
2345            bootstrap_ci_functional_logistic(&data, &y_bin, None, 2, 50, 0.05, 123, 25, 1e-6);
2346        assert!(
2347            result.is_some(),
2348            "bootstrap_ci_functional_logistic should succeed"
2349        );
2350        let ci = result.unwrap();
2351        assert_eq!(ci.lower.len(), 20);
2352        assert_eq!(ci.upper.len(), 20);
2353        assert_eq!(ci.center.len(), 20);
2354        assert!(ci.n_boot_success > 0);
2355    }
2356
2357    #[test]
2358    fn test_bootstrap_ci_logistic_ordering() {
2359        let (data, y_cont, _t) = generate_test_data(40, 20, 42);
2360        let y_median = {
2361            let mut sorted = y_cont.clone();
2362            sorted.sort_by(|a, b| a.partial_cmp(b).unwrap());
2363            sorted[sorted.len() / 2]
2364        };
2365        let y_bin: Vec<f64> = y_cont
2366            .iter()
2367            .map(|&v| if v >= y_median { 1.0 } else { 0.0 })
2368            .collect();
2369
2370        let ci = bootstrap_ci_functional_logistic(&data, &y_bin, None, 2, 100, 0.05, 42, 25, 1e-6)
2371            .unwrap();
2372        for j in 0..20 {
2373            assert!(
2374                ci.lower[j] <= ci.upper[j] + 1e-10,
2375                "lower <= upper at j={}",
2376                j
2377            );
2378        }
2379    }
2380
2381    // ----- fregre_basis_cv tests -----
2382
2383    #[test]
2384    fn test_fregre_basis_cv_returns_result() {
2385        let (data, y, t) = generate_test_data(30, 20, 42);
2386        let result = fregre_basis_cv(
2387            &data,
2388            &y,
2389            &t,
2390            5,
2391            None,
2392            7,
2393            &crate::smooth_basis::BasisType::Bspline { order: 4 },
2394        );
2395        assert!(result.is_some(), "fregre_basis_cv should succeed");
2396        let res = result.unwrap();
2397        assert!(res.optimal_lambda > 0.0);
2398        assert_eq!(res.cv_errors.len(), 20); // default 20 lambdas
2399        assert!(res.min_cv_error >= 0.0);
2400    }
2401
2402    #[test]
2403    fn test_fregre_basis_cv_custom_lambdas() {
2404        let (data, y, t) = generate_test_data(25, 15, 42);
2405        let lambdas = vec![0.001, 0.01, 0.1, 1.0, 10.0];
2406        let result = fregre_basis_cv(
2407            &data,
2408            &y,
2409            &t,
2410            5,
2411            Some(&lambdas),
2412            7,
2413            &crate::smooth_basis::BasisType::Bspline { order: 4 },
2414        );
2415        assert!(result.is_some());
2416        let res = result.unwrap();
2417        assert_eq!(res.lambda_values.len(), 5);
2418        assert!(lambdas.contains(&res.optimal_lambda));
2419    }
2420
2421    // ----- fregre_np_cv tests -----
2422
2423    #[test]
2424    fn test_fregre_np_cv_returns_result() {
2425        let (data, y, t) = generate_test_data(25, 15, 42);
2426        let result = fregre_np_cv(&data, &y, &t, 5, None, None);
2427        assert!(result.is_some(), "fregre_np_cv should succeed");
2428        let res = result.unwrap();
2429        assert!(res.optimal_h > 0.0);
2430        assert_eq!(res.cv_errors.len(), 20); // default 20 quantiles
2431        assert!(res.min_cv_error >= 0.0);
2432    }
2433
2434    #[test]
2435    fn test_fregre_np_cv_custom_h() {
2436        let (data, y, t) = generate_test_data(20, 10, 42);
2437        let h_vals = vec![0.1, 0.5, 1.0, 2.0];
2438        let result = fregre_np_cv(&data, &y, &t, 3, Some(&h_vals), None);
2439        assert!(result.is_some());
2440        let res = result.unwrap();
2441        assert_eq!(res.h_values.len(), 4);
2442        assert!(h_vals.contains(&res.optimal_h));
2443    }
2444}