Skip to main content

fdars_core/elastic_regression/
mod.rs

1//! Elastic regression models (alignment-integrated regression).
2//!
3//! These models from fdasrvf align curves during the regression fitting process,
4//! jointly optimizing alignment and regression coefficients.
5//!
6//! Key capabilities:
7//! - [`elastic_regression`] — Scalar-on-function regression with elastic alignment
8//! - [`elastic_logistic`] — Binary classification with elastic alignment
9//! - [`elastic_pcr`] — Principal component regression after elastic alignment
10//! - [`scalar_on_shape()`] — Scalar-on-shape regression with optional single-index link
11
12pub mod logistic;
13pub mod pcr;
14pub mod regression;
15pub mod scalar_on_shape;
16
17#[cfg(test)]
18mod tests;
19
20// Re-export all public items
21pub use logistic::{
22    elastic_logistic, elastic_logistic_with_config, predict_elastic_logistic, ElasticLogisticResult,
23};
24pub use pcr::{elastic_pcr, elastic_pcr_with_config, ElasticPcrResult};
25pub use regression::{
26    elastic_regression, elastic_regression_with_config, predict_elastic_regression,
27    ElasticRegressionResult,
28};
29pub use scalar_on_shape::{predict_scalar_on_shape, scalar_on_shape, ScalarOnShapeResult};
30
31use crate::alignment::reparameterize_curve;
32use crate::matrix::FdMatrix;
33
34// ─── Config Structs ─────────────────────────────────────────────────────────
35
36/// Configuration for [`elastic_regression`] and [`elastic_logistic`].
37#[derive(Debug, Clone, PartialEq)]
38pub struct ElasticConfig {
39    /// Number of basis functions for the beta coefficient (for elastic_regression).
40    pub ncomp_beta: usize,
41    /// Roughness penalty weight.
42    pub lambda: f64,
43    /// Maximum iterations for iterative alignment.
44    pub max_iter: usize,
45    /// Convergence tolerance.
46    pub tol: f64,
47}
48
49impl Default for ElasticConfig {
50    fn default() -> Self {
51        Self {
52            ncomp_beta: 10,
53            lambda: 0.0,
54            max_iter: 20,
55            tol: 1e-4,
56        }
57    }
58}
59
60/// Configuration for [`elastic_pcr`].
61#[derive(Debug, Clone, PartialEq)]
62pub struct ElasticPcrConfig {
63    /// Number of principal components to retain.
64    pub ncomp: usize,
65    /// PCA method (vertical, horizontal, or joint).
66    pub pca_method: PcaMethod,
67    /// Roughness penalty weight.
68    pub lambda: f64,
69    /// Maximum iterations for Karcher mean.
70    pub max_iter: usize,
71    /// Convergence tolerance for Karcher mean.
72    pub tol: f64,
73}
74
75impl Default for ElasticPcrConfig {
76    fn default() -> Self {
77        Self {
78            ncomp: 3,
79            pca_method: PcaMethod::Vertical,
80            lambda: 0.0,
81            max_iter: 20,
82            tol: 1e-4,
83        }
84    }
85}
86
87/// Configuration for [`scalar_on_shape()`].
88#[derive(Debug, Clone, PartialEq)]
89pub struct ScalarOnShapeConfig {
90    /// Number of Fourier basis functions for the beta representation.
91    pub nbasis: usize,
92    /// Roughness penalty weight for beta.
93    pub lambda: f64,
94    /// Penalty derivative order.
95    pub lfd_order: usize,
96    /// Index function method.
97    pub index_method: IndexMethod,
98    /// Polynomial degree for g (intercept link function).
99    pub g_degree: usize,
100    /// Maximum outer iterations (alternating beta, h, g).
101    pub max_iter_outer: usize,
102    /// Maximum inner iterations (beta estimation with alignment).
103    pub max_iter_inner: usize,
104    /// Convergence tolerance.
105    pub tol: f64,
106    /// DP alignment penalty.
107    pub dp_lambda: f64,
108}
109
110impl Default for ScalarOnShapeConfig {
111    fn default() -> Self {
112        Self {
113            nbasis: 11,
114            lambda: 1e-3,
115            lfd_order: 2,
116            index_method: IndexMethod::Identity,
117            g_degree: 2,
118            max_iter_outer: 10,
119            max_iter_inner: 15,
120            tol: 1e-4,
121            dp_lambda: 0.0,
122        }
123    }
124}
125
126// ─── Types ──────────────────────────────────────────────────────────────────
127
128/// PCA method for elastic PCR.
129#[derive(Debug, Clone, Copy, PartialEq)]
130#[non_exhaustive]
131pub enum PcaMethod {
132    Vertical,
133    Horizontal,
134    Joint,
135}
136
137/// Index function method for scalar-on-shape regression.
138///
139/// Controls the link between the shape score and the response variable.
140#[derive(Debug, Clone, PartialEq)]
141#[non_exhaustive]
142pub enum IndexMethod {
143    /// Identity: h(z) = z (standard ScoSh).
144    Identity,
145    /// Polynomial single-index: h(z) = sum of a_k z^k (SI-ScoSh).
146    Polynomial(usize),
147    /// Nadaraya-Watson kernel estimate with the given bandwidth.
148    NadarayaWatson(f64),
149}
150
151// ─── Shared Helpers ────────────────────────────────────────────────────────
152
153/// Apply warping functions to SRSFs, producing aligned SRSFs with sqrt(γ') factor.
154pub(super) fn apply_warps_to_srsfs(
155    q_all: &FdMatrix,
156    gammas: &FdMatrix,
157    argvals: &[f64],
158) -> FdMatrix {
159    let (n, m) = q_all.shape();
160    let h = (argvals[m - 1] - argvals[0]) / (m - 1) as f64;
161    let mut q_aligned = FdMatrix::zeros(n, m);
162    for i in 0..n {
163        let qi: Vec<f64> = (0..m).map(|j| q_all[(i, j)]).collect();
164        let gam: Vec<f64> = (0..m).map(|j| gammas[(i, j)]).collect();
165        let q_warped = reparameterize_curve(&qi, argvals, &gam);
166        let gam_deriv = crate::helpers::gradient_uniform(&gam, h);
167        for j in 0..m {
168            q_aligned[(i, j)] = q_warped[j] * gam_deriv[j].max(0.0).sqrt();
169        }
170    }
171    q_aligned
172}
173
174/// Initialize warping functions to identity (γ_i(t) = t).
175pub(super) fn init_identity_warps(n: usize, argvals: &[f64]) -> FdMatrix {
176    let m = argvals.len();
177    let mut gammas = FdMatrix::zeros(n, m);
178    for i in 0..n {
179        for j in 0..m {
180            gammas[(i, j)] = argvals[j];
181        }
182    }
183    gammas
184}
185
186/// Compute fitted values: ŷ_i = α + ∫ q_aligned_i · β · w dt.
187pub(super) fn srsf_fitted_values(
188    q_aligned: &FdMatrix,
189    beta: &[f64],
190    weights: &[f64],
191    alpha: f64,
192) -> Vec<f64> {
193    let (n, m) = q_aligned.shape();
194    let mut fitted = vec![0.0; n];
195    for i in 0..n {
196        fitted[i] = alpha;
197        for j in 0..m {
198            fitted[i] += q_aligned[(i, j)] * beta[j] * weights[j];
199        }
200    }
201    fitted
202}
203
204/// Check relative convergence of β.
205pub(super) fn beta_converged(beta_new: &[f64], beta_old: &[f64], tol: f64) -> bool {
206    let diff: f64 = beta_new
207        .iter()
208        .zip(beta_old.iter())
209        .map(|(&a, &b)| (a - b).powi(2))
210        .sum::<f64>()
211        .sqrt();
212    let norm: f64 = beta_old
213        .iter()
214        .map(|&b| b * b)
215        .sum::<f64>()
216        .sqrt()
217        .max(1e-10);
218    diff / norm < tol
219}