Skip to main content

kriging_rs/
projected.rs

1//! Projected (planar) coordinates with Euclidean distance and 2-D anisotropy.
2//!
3//! The rest of the crate is geographic (`GeoCoord` + Haversine). This module adds:
4//!
5//! - [`ProjectedCoord`] — a Cartesian `(x, y)` pair in any linear unit (meters, km, etc.).
6//! - [`euclidean_distance`] and [`euclidean_distance_squared`].
7//! - [`Anisotropy2D`] — an angle + ratio that deforms distance to model direction-dependent
8//!   range (e.g. stronger correlation along a 30° azimuth than perpendicular to it).
9//!
10//! A companion [`ProjectedKrigingModel`] performs ordinary kriging on projected coordinates,
11//! optionally with an anisotropy. Geographic data can be converted to projected data with
12//! [`ProjectedCoord::equirectangular`] for small-area work.
13
14use nalgebra::{DMatrix, DVector, Dyn, linalg::LU};
15#[cfg(not(target_arch = "wasm32"))]
16use rayon::prelude::*;
17
18use std::num::NonZeroUsize;
19
20use crate::Real;
21use crate::distance::GeoCoord;
22use crate::error::KrigingError;
23use crate::kriging::binomial::{
24    BINOMIAL_CALIBRATION_VERSION, BinomialBuildNotes, BinomialCalibratedResult, BinomialPrediction,
25    BinomialPrior, HeteroskedasticBinomialConfig, logit_observation_variance_empirical_bayes,
26};
27
28use crate::kriging::ordinary::{Prediction, kriging_diagonal_jitter};
29use crate::utils::{Probability, logistic, logit};
30use crate::variogram::empirical::{EmpiricalVariogram, PositiveReal};
31use crate::variogram::models::VariogramModel;
32
33/// A planar (projected) coordinate. Units are whatever the user supplies, as long as they
34/// match the variogram's range units.
35#[derive(Debug, Clone, Copy, PartialEq)]
36pub struct ProjectedCoord {
37    pub x: Real,
38    pub y: Real,
39}
40
41impl ProjectedCoord {
42    pub fn new(x: Real, y: Real) -> Self {
43        Self { x, y }
44    }
45
46    /// Equirectangular projection of a [`GeoCoord`] around a reference point. Returns `(x, y)`
47    /// in kilometers. Accurate for small areas; not suitable for global-scale data.
48    ///
49    /// `x = R · (lon − lon₀) · cos(lat₀)`, `y = R · (lat − lat₀)` with `R = 6371 km`.
50    pub fn equirectangular(coord: GeoCoord, reference: GeoCoord) -> Self {
51        const EARTH_R_KM: Real = 6_371.0;
52        let cos_lat0 = reference.lat().to_radians().cos();
53        let dlat = (coord.lat() - reference.lat()).to_radians();
54        let dlon = (coord.lon() - reference.lon()).to_radians();
55        Self {
56            x: EARTH_R_KM * dlon * cos_lat0,
57            y: EARTH_R_KM * dlat,
58        }
59    }
60}
61
62/// Squared Euclidean distance. Cheaper than [`euclidean_distance`] when only ordering matters
63/// or when the caller will square-root afterwards.
64#[inline]
65pub fn euclidean_distance_squared(a: ProjectedCoord, b: ProjectedCoord) -> Real {
66    let dx = a.x - b.x;
67    let dy = a.y - b.y;
68    dx * dx + dy * dy
69}
70
71/// Straight-line distance between two planar points.
72#[inline]
73pub fn euclidean_distance(a: ProjectedCoord, b: ProjectedCoord) -> Real {
74    euclidean_distance_squared(a, b).sqrt()
75}
76
77/// Geometric anisotropy in 2-D: direction-dependent range.
78///
79/// `major_angle_deg` is the orientation (counter-clockwise from the +x axis) of the axis of
80/// **maximum** continuity, in degrees. `range_ratio = minor_range / major_range ∈ (0, 1]` is
81/// the ratio of minor-axis to major-axis correlation range. `range_ratio = 1` is isotropic.
82///
83/// The effective distance between two points in anisotropy space is
84///
85/// ```text
86///   h' = sqrt((h_major)² + (h_minor / range_ratio)²)
87/// ```
88///
89/// so the minor-axis separation is scaled up (shortening correlation) by `1 / range_ratio`.
90#[derive(Debug, Clone, Copy, PartialEq)]
91pub struct Anisotropy2D {
92    pub major_angle_deg: Real,
93    pub range_ratio: Real,
94}
95
96impl Anisotropy2D {
97    /// Construct a new anisotropy. Returns an error if `range_ratio` is not in `(0, 1]`
98    /// or is non-finite.
99    pub fn new(major_angle_deg: Real, range_ratio: Real) -> Result<Self, KrigingError> {
100        if !range_ratio.is_finite() || range_ratio <= 0.0 || range_ratio > 1.0 {
101            return Err(KrigingError::InvalidInput(format!(
102                "range_ratio must be in (0, 1], got {range_ratio}"
103            )));
104        }
105        if !major_angle_deg.is_finite() {
106            return Err(KrigingError::InvalidInput(
107                "major_angle_deg must be finite".to_string(),
108            ));
109        }
110        Ok(Self {
111            major_angle_deg,
112            range_ratio,
113        })
114    }
115
116    /// Isotropic anisotropy (no deformation). Equivalent to unmodified Euclidean distance.
117    pub fn isotropic() -> Self {
118        Self {
119            major_angle_deg: 0.0,
120            range_ratio: 1.0,
121        }
122    }
123
124    /// Effective (anisotropy-adjusted) distance between two planar points.
125    pub fn distance(self, a: ProjectedCoord, b: ProjectedCoord) -> Real {
126        let dx = b.x - a.x;
127        let dy = b.y - a.y;
128        let theta = self.major_angle_deg.to_radians();
129        let cos_t = theta.cos();
130        let sin_t = theta.sin();
131        // Rotate so the major axis aligns with +x.
132        let h_major = dx * cos_t + dy * sin_t;
133        let h_minor = -dx * sin_t + dy * cos_t;
134        let inv_r = 1.0 / self.range_ratio;
135        (h_major * h_major + (h_minor * inv_r) * (h_minor * inv_r)).sqrt()
136    }
137}
138
139/// Directional empirical variogram configuration for planar data.
140///
141/// Pairs are accepted only when the vector from `i → j` lies within `tolerance_deg` of
142/// `azimuth_deg` (or its 180°-mirror — the variogram is symmetric under direction reversal).
143/// `azimuth_deg` is measured counter-clockwise from +x (i.e. "east"). Use
144/// `azimuth_deg = 0` for east–west, `90` for north–south.
145#[derive(Debug, Clone, Copy, PartialEq)]
146pub struct DirectionalConfig {
147    pub azimuth_deg: Real,
148    pub tolerance_deg: Real,
149}
150
151impl DirectionalConfig {
152    pub fn new(azimuth_deg: Real, tolerance_deg: Real) -> Result<Self, KrigingError> {
153        if !tolerance_deg.is_finite() || tolerance_deg <= 0.0 || tolerance_deg > 90.0 {
154            return Err(KrigingError::InvalidInput(format!(
155                "tolerance_deg must be in (0, 90], got {tolerance_deg}"
156            )));
157        }
158        if !azimuth_deg.is_finite() {
159            return Err(KrigingError::InvalidInput(
160                "azimuth_deg must be finite".to_string(),
161            ));
162        }
163        Ok(Self {
164            azimuth_deg,
165            tolerance_deg,
166        })
167    }
168}
169
170/// Compute a directional empirical variogram on projected coordinates. Only pairs whose
171/// separation vector falls within the angular tolerance around `direction.azimuth_deg`
172/// contribute. See [`compute_empirical_variogram`](crate::compute_empirical_variogram) for
173/// the omnidirectional, geographic counterpart.
174pub fn compute_directional_empirical_variogram(
175    coords: &[ProjectedCoord],
176    values: &[Real],
177    max_distance: PositiveReal,
178    n_bins: NonZeroUsize,
179    direction: DirectionalConfig,
180) -> Result<EmpiricalVariogram, KrigingError> {
181    if coords.len() != values.len() {
182        return Err(KrigingError::DimensionMismatch(format!(
183            "coords ({}) and values ({}) must have equal length",
184            coords.len(),
185            values.len()
186        )));
187    }
188    if coords.is_empty() {
189        return Err(KrigingError::InsufficientData(1));
190    }
191    let n_bins = n_bins.get();
192    let bin_width = max_distance.get() / n_bins as Real;
193    let mut dist_sums = vec![0.0 as Real; n_bins];
194    let mut semi_sums = vec![0.0 as Real; n_bins];
195    let mut counts = vec![0usize; n_bins];
196
197    let target = direction.azimuth_deg.to_radians();
198    let tol = direction.tolerance_deg.to_radians();
199
200    for i in 0..coords.len() {
201        for j in (i + 1)..coords.len() {
202            let dx = coords[j].x - coords[i].x;
203            let dy = coords[j].y - coords[i].y;
204            let d = (dx * dx + dy * dy).sqrt();
205            if d <= 0.0 || d > max_distance.get() {
206                continue;
207            }
208            let angle = dy.atan2(dx);
209            // Fold into [-π/2, π/2] since the variogram is symmetric under direction reversal.
210            let mut diff = angle - target;
211            while diff > std::f64::consts::FRAC_PI_2 as Real {
212                diff -= std::f64::consts::PI as Real;
213            }
214            while diff < -(std::f64::consts::FRAC_PI_2 as Real) {
215                diff += std::f64::consts::PI as Real;
216            }
217            if diff.abs() > tol {
218                continue;
219            }
220            let g = 0.5 * (values[i] - values[j]).powi(2);
221            let mut bin = (d / bin_width).floor() as usize;
222            if bin >= n_bins {
223                bin = n_bins - 1;
224            }
225            dist_sums[bin] += d;
226            semi_sums[bin] += g;
227            counts[bin] += 1;
228        }
229    }
230
231    let mut distances = Vec::new();
232    let mut semivariances = Vec::new();
233    let mut n_pairs = Vec::new();
234    for i in 0..n_bins {
235        if counts[i] == 0 {
236            continue;
237        }
238        distances.push(dist_sums[i] / counts[i] as Real);
239        semivariances.push(semi_sums[i] / counts[i] as Real);
240        n_pairs.push(counts[i]);
241    }
242    if distances.is_empty() {
243        return Err(KrigingError::FittingError(
244            "no pairs in selected direction/distance range".to_string(),
245        ));
246    }
247    Ok(EmpiricalVariogram {
248        distances,
249        semivariances,
250        n_pairs,
251    })
252}
253
254/// Dataset of projected coordinates and paired values. Construction validates lengths.
255#[derive(Debug, Clone)]
256pub struct ProjectedDataset {
257    pub coords: Vec<ProjectedCoord>,
258    pub values: Vec<Real>,
259}
260
261impl ProjectedDataset {
262    pub fn new(coords: Vec<ProjectedCoord>, values: Vec<Real>) -> Result<Self, KrigingError> {
263        if coords.len() != values.len() {
264            return Err(KrigingError::DimensionMismatch(format!(
265                "coords ({}) and values ({}) must have equal length",
266                coords.len(),
267                values.len()
268            )));
269        }
270        if coords.is_empty() {
271            return Err(KrigingError::InsufficientData(1));
272        }
273        Ok(Self { coords, values })
274    }
275}
276
277/// Ordinary kriging on projected (planar) coordinates with optional 2-D geometric
278/// anisotropy. The variogram's range must be expressed in the same linear units as the
279/// coordinates (e.g. km, m).
280#[derive(Debug)]
281pub struct ProjectedKrigingModel {
282    coords: Vec<ProjectedCoord>,
283    values: Vec<Real>,
284    variogram: VariogramModel,
285    anisotropy: Anisotropy2D,
286    cov_at_zero: Real,
287    system: DMatrix<Real>,
288    system_lu: LU<Real, Dyn, Dyn>,
289}
290
291impl Clone for ProjectedKrigingModel {
292    fn clone(&self) -> Self {
293        let system = self.system.clone();
294        let system_lu = system.clone().lu();
295        Self {
296            coords: self.coords.clone(),
297            values: self.values.clone(),
298            variogram: self.variogram,
299            anisotropy: self.anisotropy,
300            cov_at_zero: self.cov_at_zero,
301            system,
302            system_lu,
303        }
304    }
305}
306
307impl ProjectedKrigingModel {
308    /// Build a projected ordinary kriging model. Pass [`Anisotropy2D::isotropic`] to disable
309    /// anisotropy.
310    pub fn new(
311        dataset: ProjectedDataset,
312        variogram: VariogramModel,
313        anisotropy: Anisotropy2D,
314    ) -> Result<Self, KrigingError> {
315        Self::new_with_extra_diagonal_internal(dataset, variogram, anisotropy, &[])
316    }
317
318    /// Like [`new`](Self::new) but adds `extra` (length `n`) to each main-diagonal covariance
319    /// term, modeling observation-specific (non-spatial) noise in addition to jitter.
320    pub fn new_with_extra_diagonal(
321        dataset: ProjectedDataset,
322        variogram: VariogramModel,
323        anisotropy: Anisotropy2D,
324        extra: Vec<Real>,
325    ) -> Result<Self, KrigingError> {
326        if !extra.is_empty() && extra.len() != dataset.coords.len() {
327            return Err(KrigingError::InvalidInput(
328                "extra observation diagonal must be empty (homoscedastic) or the same length as the dataset"
329                    .to_string(),
330            ));
331        }
332        for &v in &extra {
333            if !v.is_finite() || v < 0.0 {
334                return Err(KrigingError::InvalidInput(
335                    "observation diagonal entries must be finite and non-negative".to_string(),
336                ));
337            }
338        }
339        Self::new_with_extra_diagonal_internal(dataset, variogram, anisotropy, &extra)
340    }
341
342    fn new_with_extra_diagonal_internal(
343        dataset: ProjectedDataset,
344        variogram: VariogramModel,
345        anisotropy: Anisotropy2D,
346        extra: &[Real],
347    ) -> Result<Self, KrigingError> {
348        let coords = dataset.coords;
349        let values = dataset.values;
350        let n = coords.len();
351        if !extra.is_empty() && extra.len() != n {
352            return Err(KrigingError::InvalidInput(
353                "internal: extra length mismatch for projected kriging".to_string(),
354            ));
355        }
356        let mut system = DMatrix::from_element(n + 1, n + 1, 0.0);
357        let diag_eps = kriging_diagonal_jitter(n, variogram);
358        for i in 0..n {
359            for j in i..n {
360                let d = anisotropy.distance(coords[i], coords[j]);
361                let mut cov = variogram.covariance(d);
362                if i == j {
363                    cov += diag_eps;
364                    if let Some(&dextra) = extra.get(i) {
365                        cov += dextra;
366                    }
367                }
368                system[(i, j)] = cov;
369                system[(j, i)] = cov;
370            }
371            system[(i, n)] = 1.0;
372            system[(n, i)] = 1.0;
373        }
374        system[(n, n)] = 0.0;
375
376        let system_lu = system.clone().lu();
377        let mut probe = DVector::from_element(n + 1, 0.0);
378        probe[n] = 1.0;
379        if system_lu.solve(&probe).is_none() {
380            return Err(KrigingError::MatrixError(
381                "could not factorize projected kriging system".to_string(),
382            ));
383        }
384
385        Ok(Self {
386            coords,
387            values,
388            variogram,
389            anisotropy,
390            cov_at_zero: variogram.covariance(0.0),
391            system,
392            system_lu,
393        })
394    }
395
396    pub fn anisotropy(&self) -> Anisotropy2D {
397        self.anisotropy
398    }
399
400    pub fn predict(&self, coord: ProjectedCoord) -> Result<Prediction, KrigingError> {
401        let mut rhs = DVector::from_element(self.coords.len() + 1, 0.0);
402        self.predict_with_rhs(coord, &mut rhs)
403    }
404
405    pub fn predict_batch(
406        &self,
407        coords: &[ProjectedCoord],
408    ) -> Result<Vec<Prediction>, KrigingError> {
409        #[cfg(not(target_arch = "wasm32"))]
410        {
411            let n = self.coords.len();
412            coords
413                .par_iter()
414                .map(|c| {
415                    let mut rhs = DVector::from_element(n + 1, 0.0);
416                    self.predict_with_rhs(*c, &mut rhs)
417                })
418                .collect()
419        }
420        #[cfg(target_arch = "wasm32")]
421        {
422            let mut rhs = DVector::from_element(self.coords.len() + 1, 0.0);
423            let mut out = Vec::with_capacity(coords.len());
424            for &c in coords {
425                out.push(self.predict_with_rhs(c, &mut rhs)?);
426            }
427            Ok(out)
428        }
429    }
430
431    fn predict_with_rhs(
432        &self,
433        coord: ProjectedCoord,
434        rhs: &mut DVector<Real>,
435    ) -> Result<Prediction, KrigingError> {
436        let n = self.coords.len();
437        for i in 0..n {
438            let d = self.anisotropy.distance(self.coords[i], coord);
439            rhs[i] = self.variogram.covariance(d);
440        }
441        rhs[n] = 1.0;
442        let sol = self.system_lu.solve(rhs).ok_or_else(|| {
443            KrigingError::MatrixError("could not solve projected kriging system".to_string())
444        })?;
445        let mut value: Real = 0.0;
446        let mut cov_dot: Real = 0.0;
447        for i in 0..n {
448            value += sol[i] * self.values[i];
449            cov_dot += sol[i] * rhs[i];
450        }
451        let mu = sol[n];
452        let variance = (self.cov_at_zero - cov_dot - mu).max(0.0);
453        Ok(Prediction { value, variance })
454    }
455}
456
457// ---------------------------------------------------------------------------
458// Binomial kriging on projected coordinates
459// ---------------------------------------------------------------------------
460
461/// A single binomial observation on projected (planar) coordinates.
462///
463/// This is the [`crate::BinomialObservation`] equivalent for
464/// [`BinomialProjectedKrigingModel`] / planar prevalence mapping.
465#[derive(Debug, Clone, Copy)]
466pub struct ProjectedBinomialObservation {
467    coord: ProjectedCoord,
468    successes: u32,
469    trials: u32,
470}
471
472impl ProjectedBinomialObservation {
473    /// Validates `trials > 0` and `successes <= trials`.
474    pub fn new(coord: ProjectedCoord, successes: u32, trials: u32) -> Result<Self, KrigingError> {
475        if trials == 0 {
476            return Err(KrigingError::InvalidBinomialData(
477                "trials must be greater than 0".to_string(),
478            ));
479        }
480        if successes > trials {
481            return Err(KrigingError::InvalidBinomialData(format!(
482                "successes ({}) cannot exceed trials ({})",
483                successes, trials
484            )));
485        }
486        Ok(Self {
487            coord,
488            successes,
489            trials,
490        })
491    }
492
493    #[inline]
494    pub fn coord(self) -> ProjectedCoord {
495        self.coord
496    }
497
498    #[inline]
499    pub fn successes(self) -> u32 {
500        self.successes
501    }
502
503    #[inline]
504    pub fn trials(self) -> u32 {
505        self.trials
506    }
507
508    pub fn smoothed_probability_with_prior(&self, prior: BinomialPrior) -> Real {
509        let s = self.successes as Real;
510        let n = self.trials as Real;
511        (s + prior.alpha()) / (n + prior.alpha() + prior.beta())
512    }
513
514    pub fn smoothed_logit_with_prior(&self, prior: BinomialPrior) -> Real {
515        let p = self.smoothed_probability_with_prior(prior);
516        logit(Probability::from_known_in_range(p))
517    }
518}
519
520#[inline]
521fn delta_prevalence_variance(prevalence: Real, logit_variance: Real) -> Real {
522    let factor = prevalence * (1.0 - prevalence);
523    factor * factor * logit_variance.max(0.0)
524}
525
526fn binomial_prediction_from(pred: Prediction) -> BinomialPrediction {
527    let prevalence = logistic(pred.value);
528    BinomialPrediction {
529        prevalence,
530        logit_value: pred.value,
531        variance: pred.variance,
532        prevalence_variance: delta_prevalence_variance(prevalence, pred.variance),
533    }
534}
535
536/// Binomial kriging on projected (planar) coordinates with optional 2-D
537/// geometric anisotropy.
538///
539/// This is the planar/anisotropic counterpart of
540/// [`crate::BinomialKrigingModel`]. Kriging happens on the **logit** scale, with
541/// per-site logit observation variance (calibrated binomial path) when built from
542/// counts. See [`BinomialBuildNotes`].
543///
544/// Use [`Anisotropy2D::isotropic`] to disable anisotropy.
545#[derive(Debug, Clone)]
546pub struct BinomialProjectedKrigingModel {
547    inner: ProjectedKrigingModel,
548}
549
550/// Planar / anisotropic binomial fit: model + [`BinomialBuildNotes`].
551pub type ProjectedBinomialFit = BinomialCalibratedResult<BinomialProjectedKrigingModel>;
552
553impl BinomialProjectedKrigingModel {
554    /// Build a binomial projected kriging model from `(coord, successes, trials)` with the
555    /// default `Beta(1, 1)` prior.
556    pub fn new(
557        observations: Vec<ProjectedBinomialObservation>,
558        variogram: VariogramModel,
559        anisotropy: Anisotropy2D,
560    ) -> Result<ProjectedBinomialFit, KrigingError> {
561        Self::new_with_prior(
562            observations,
563            variogram,
564            anisotropy,
565            BinomialPrior::default(),
566        )
567    }
568
569    /// As [`new`](Self::new), with an explicit Beta prior. Uses heteroskedastic (calibrated) conditioning.
570    pub fn new_with_prior(
571        observations: Vec<ProjectedBinomialObservation>,
572        variogram: VariogramModel,
573        anisotropy: Anisotropy2D,
574        prior: BinomialPrior,
575    ) -> Result<ProjectedBinomialFit, KrigingError> {
576        if observations.len() < 2 {
577            return Err(KrigingError::InsufficientData(2));
578        }
579        let coords: Vec<ProjectedCoord> = observations.iter().map(|o| o.coord()).collect();
580        let logits: Vec<Real> = observations
581            .iter()
582            .map(|o| o.smoothed_logit_with_prior(prior))
583            .collect();
584        let base: Vec<Real> = observations
585            .iter()
586            .map(|o| logit_observation_variance_empirical_bayes(prior, o.successes(), o.trials()))
587            .map(|v| v.max(HeteroskedasticBinomialConfig::default().min_logit_observation_variance))
588            .collect();
589        let n_tries = HeteroskedasticBinomialConfig::default()
590            .max_build_attempts
591            .max(1);
592        let config = HeteroskedasticBinomialConfig::default();
593        let mut last_err: Option<KrigingError> = None;
594        let mut inflation = 1.0 as Real;
595        for attempt in 0..n_tries {
596            let extra: Vec<Real> = base
597                .iter()
598                .map(|&v| (v * inflation).max(config.min_logit_observation_variance))
599                .collect();
600            let dataset = ProjectedDataset::new(coords.clone(), logits.clone())?;
601            match ProjectedKrigingModel::new_with_extra_diagonal(
602                dataset, variogram, anisotropy, extra,
603            ) {
604                Ok(inner) => {
605                    return Ok(BinomialCalibratedResult {
606                        model: Self { inner },
607                        notes: BinomialBuildNotes {
608                            calibration_version: BINOMIAL_CALIBRATION_VERSION,
609                            logit_inflation: inflation,
610                            n_build_attempts: attempt + 1,
611                            prior,
612                            zero_trial_dropped_indices: Vec::new(),
613                            from_precomputed_logits_only: false,
614                        },
615                    });
616                }
617                Err(e) => {
618                    last_err = Some(e);
619                }
620            }
621            inflation *= 2.0 as Real;
622        }
623        Err(last_err.unwrap_or_else(|| {
624            KrigingError::MatrixError("binomial projected kriging build failed".to_string())
625        }))
626    }
627
628    /// Pre-computed logits with no per-trial data (no observation variance on the diagonal
629    /// except variogram and jitter). See [`crate::BinomialKrigingModel::from_precomputed_logits`].
630    pub fn from_precomputed_logits(
631        coords: Vec<ProjectedCoord>,
632        logits: Vec<Real>,
633        variogram: VariogramModel,
634        anisotropy: Anisotropy2D,
635    ) -> Result<ProjectedBinomialFit, KrigingError> {
636        if logits.iter().any(|v| !v.is_finite()) {
637            return Err(KrigingError::InvalidInput(
638                "logits must all be finite (no NaN/inf)".to_string(),
639            ));
640        }
641        let dataset = ProjectedDataset::new(coords, logits)?;
642        let inner = ProjectedKrigingModel::new(dataset, variogram, anisotropy)?;
643        Ok(BinomialCalibratedResult {
644            model: Self { inner },
645            notes: BinomialBuildNotes {
646                calibration_version: BINOMIAL_CALIBRATION_VERSION,
647                logit_inflation: 1.0,
648                n_build_attempts: 1,
649                prior: BinomialPrior::default(),
650                zero_trial_dropped_indices: Vec::new(),
651                from_precomputed_logits_only: true,
652            },
653        })
654    }
655
656    /// Pre-computed logit field with per-site logit observation variances and stability policy.
657    pub fn from_precomputed_logits_with_logit_observation_variances(
658        coords: Vec<ProjectedCoord>,
659        logits: Vec<Real>,
660        variogram: VariogramModel,
661        anisotropy: Anisotropy2D,
662        base_logit_observation_variance: Vec<Real>,
663        config: HeteroskedasticBinomialConfig,
664        prior_for_notes: BinomialPrior,
665    ) -> Result<ProjectedBinomialFit, KrigingError> {
666        if logits.len() != coords.len() {
667            return Err(KrigingError::DimensionMismatch(
668                "coords and logits must have equal length".to_string(),
669            ));
670        }
671        if base_logit_observation_variance.len() != coords.len() {
672            return Err(KrigingError::InvalidInput(
673                "logit observation variance must match coords length".to_string(),
674            ));
675        }
676        if logits.iter().any(|v| !v.is_finite()) {
677            return Err(KrigingError::InvalidInput(
678                "logits must all be finite (no NaN/inf)".to_string(),
679            ));
680        }
681        for &v in &base_logit_observation_variance {
682            if !v.is_finite() || v < 0.0 {
683                return Err(KrigingError::InvalidInput(
684                    "logit observation variances must be finite and non-negative".to_string(),
685                ));
686            }
687        }
688        let n_tries = config.max_build_attempts.max(1);
689        let mut last_err: Option<KrigingError> = None;
690        let mut inflation = 1.0 as Real;
691        for attempt in 0..n_tries {
692            let extra: Vec<Real> = base_logit_observation_variance
693                .iter()
694                .map(|&v| (v * inflation).max(config.min_logit_observation_variance))
695                .collect();
696            let dataset = ProjectedDataset::new(coords.clone(), logits.clone())?;
697            match ProjectedKrigingModel::new_with_extra_diagonal(
698                dataset, variogram, anisotropy, extra,
699            ) {
700                Ok(inner) => {
701                    return Ok(BinomialCalibratedResult {
702                        model: Self { inner },
703                        notes: BinomialBuildNotes {
704                            calibration_version: BINOMIAL_CALIBRATION_VERSION,
705                            logit_inflation: inflation,
706                            n_build_attempts: attempt + 1,
707                            prior: prior_for_notes,
708                            zero_trial_dropped_indices: Vec::new(),
709                            from_precomputed_logits_only: false,
710                        },
711                    });
712                }
713                Err(e) => {
714                    last_err = Some(e);
715                }
716            }
717            inflation *= 2.0 as Real;
718        }
719        Err(last_err.unwrap_or_else(|| {
720            KrigingError::MatrixError(
721                "from_precomputed (projected) with observation variances: build failed".to_string(),
722            )
723        }))
724    }
725
726    pub fn anisotropy(&self) -> Anisotropy2D {
727        self.inner.anisotropy()
728    }
729
730    pub fn predict(&self, coord: ProjectedCoord) -> Result<BinomialPrediction, KrigingError> {
731        let pred = self.inner.predict(coord)?;
732        Ok(binomial_prediction_from(pred))
733    }
734
735    pub fn predict_batch(
736        &self,
737        coords: &[ProjectedCoord],
738    ) -> Result<Vec<BinomialPrediction>, KrigingError> {
739        let preds = self.inner.predict_batch(coords)?;
740        Ok(preds.into_iter().map(binomial_prediction_from).collect())
741    }
742}
743
744#[cfg(test)]
745mod tests {
746    use super::*;
747    use crate::variogram::models::VariogramType;
748
749    #[test]
750    fn euclidean_distance_matches_pythagoras() {
751        let a = ProjectedCoord::new(0.0, 0.0);
752        let b = ProjectedCoord::new(3.0, 4.0);
753        assert!((euclidean_distance(a, b) - 5.0).abs() < 1e-6);
754    }
755
756    #[test]
757    fn anisotropy_rejects_invalid_ratio() {
758        assert!(Anisotropy2D::new(0.0, 0.0).is_err());
759        assert!(Anisotropy2D::new(0.0, -1.0).is_err());
760        assert!(Anisotropy2D::new(0.0, 1.5).is_err());
761        assert!(Anisotropy2D::new(Real::NAN, 0.5).is_err());
762    }
763
764    #[test]
765    fn isotropic_distance_equals_euclidean() {
766        let iso = Anisotropy2D::isotropic();
767        let a = ProjectedCoord::new(1.0, 2.0);
768        let b = ProjectedCoord::new(4.0, 6.0);
769        let want = euclidean_distance(a, b);
770        assert!((iso.distance(a, b) - want).abs() < 1e-6);
771    }
772
773    #[test]
774    fn anisotropy_stretches_minor_axis_distance() {
775        // Major axis along +x (angle 0), ratio 0.5 => distances along +y are doubled.
776        let aniso = Anisotropy2D::new(0.0, 0.5).unwrap();
777        let origin = ProjectedCoord::new(0.0, 0.0);
778        let along_x = ProjectedCoord::new(1.0, 0.0);
779        let along_y = ProjectedCoord::new(0.0, 1.0);
780        assert!((aniso.distance(origin, along_x) - 1.0).abs() < 1e-6);
781        assert!((aniso.distance(origin, along_y) - 2.0).abs() < 1e-6);
782    }
783
784    #[test]
785    fn anisotropy_rotation_invariant_at_45_with_isotropic_ratio() {
786        let aniso = Anisotropy2D::new(37.0, 1.0).unwrap();
787        let a = ProjectedCoord::new(0.0, 0.0);
788        let b = ProjectedCoord::new(3.0, 4.0);
789        assert!((aniso.distance(a, b) - 5.0).abs() < 1e-5);
790    }
791
792    #[test]
793    fn projected_model_recovers_training_value_at_collocated_point() {
794        let coords = vec![
795            ProjectedCoord::new(0.0, 0.0),
796            ProjectedCoord::new(0.0, 10.0),
797            ProjectedCoord::new(10.0, 0.0),
798        ];
799        let values = vec![10.0, 20.0, 15.0];
800        let variogram = VariogramModel::new(0.01, 5.0, 30.0, VariogramType::Exponential).unwrap();
801        let dataset = ProjectedDataset::new(coords.clone(), values).unwrap();
802        let model = ProjectedKrigingModel::new(dataset, variogram, Anisotropy2D::isotropic())
803            .expect("model");
804        let pred = model.predict(coords[0]).expect("prediction");
805        assert!((pred.value - 10.0).abs() < 1e-3);
806        assert!(pred.variance >= 0.0);
807    }
808
809    #[test]
810    fn anisotropy_makes_along_axis_prediction_closer() {
811        // Two stations at (±5, 0); major axis along x with small ratio makes the x-neighbors
812        // more influential than the same-distance y-direction would be.
813        let coords = vec![
814            ProjectedCoord::new(-5.0, 0.0),
815            ProjectedCoord::new(5.0, 0.0),
816            ProjectedCoord::new(0.0, 5.0),
817            ProjectedCoord::new(0.0, -5.0),
818        ];
819        let values = vec![10.0, 10.0, 30.0, 30.0];
820        let variogram = VariogramModel::new(0.01, 1.0, 20.0, VariogramType::Exponential).unwrap();
821        let dataset = ProjectedDataset::new(coords.clone(), values.clone()).unwrap();
822        // Strong anisotropy: along-x continuity much larger than along-y.
823        let iso_model = ProjectedKrigingModel::new(
824            ProjectedDataset::new(coords.clone(), values.clone()).unwrap(),
825            variogram,
826            Anisotropy2D::isotropic(),
827        )
828        .expect("iso");
829        let aniso = Anisotropy2D::new(0.0, 0.2).unwrap();
830        let aniso_model = ProjectedKrigingModel::new(dataset, variogram, aniso).expect("aniso");
831
832        // At origin, iso averages all stations (value ~20) while aniso overweights x-axis
833        // neighbors (value ~10).
834        let target = ProjectedCoord::new(0.0, 0.0);
835        let iso_pred = iso_model.predict(target).expect("iso pred");
836        let aniso_pred = aniso_model.predict(target).expect("aniso pred");
837        assert!(
838            aniso_pred.value < iso_pred.value,
839            "expected aniso.value={} < iso.value={}",
840            aniso_pred.value,
841            iso_pred.value
842        );
843    }
844
845    #[test]
846    fn equirectangular_small_area_matches_haversine_within_a_percent() {
847        // Project two nearby points and compare Euclidean distance to Haversine.
848        let reference = GeoCoord::try_new(10.0, 20.0).unwrap();
849        let a = GeoCoord::try_new(10.1, 20.05).unwrap();
850        let b = GeoCoord::try_new(10.2, 20.10).unwrap();
851        let pa = ProjectedCoord::equirectangular(a, reference);
852        let pb = ProjectedCoord::equirectangular(b, reference);
853        let planar = euclidean_distance(pa, pb);
854        let geo = crate::distance::haversine_distance(a, b);
855        let rel_err = (planar - geo).abs() / geo.max(1e-6);
856        assert!(rel_err < 0.01, "rel_err={rel_err}");
857    }
858
859    #[test]
860    fn directional_variogram_only_counts_pairs_in_specified_direction() {
861        // Three points arranged east-west and three arranged north-south, all distinct.
862        // An east-west direction window should only capture the east-west pairs.
863        let coords = vec![
864            ProjectedCoord { x: 0.0, y: 0.0 },
865            ProjectedCoord { x: 1.0, y: 0.0 },
866            ProjectedCoord { x: 2.0, y: 0.0 },
867            ProjectedCoord { x: 0.0, y: 10.0 },
868            ProjectedCoord { x: 0.0, y: 20.0 },
869        ];
870        let values = vec![1.0, 2.0, 3.0, 10.0, 20.0];
871        let cfg = DirectionalConfig::new(0.0, 5.0).unwrap(); // east, ±5°
872        let out = compute_directional_empirical_variogram(
873            &coords,
874            &values,
875            PositiveReal::try_new(5.0).unwrap(),
876            NonZeroUsize::new(5).unwrap(),
877            cfg,
878        )
879        .expect("directional variogram");
880        // Only the three east-west pairs should be counted (0-1, 1-2, 0-2).
881        let total_pairs: usize = out.n_pairs.iter().sum();
882        assert_eq!(total_pairs, 3);
883    }
884
885    #[test]
886    fn directional_config_rejects_invalid_tolerance() {
887        assert!(DirectionalConfig::new(0.0, 0.0).is_err());
888        assert!(DirectionalConfig::new(0.0, 91.0).is_err());
889        assert!(DirectionalConfig::new(Real::NAN, 5.0).is_err());
890    }
891
892    #[test]
893    fn binomial_projected_smoothed_logit_near_collocated_with_calibrated_diagonal() {
894        let obs = vec![
895            ProjectedBinomialObservation::new(ProjectedCoord::new(0.0, 0.0), 3, 10).unwrap(),
896            ProjectedBinomialObservation::new(ProjectedCoord::new(10.0, 0.0), 7, 10).unwrap(),
897            ProjectedBinomialObservation::new(ProjectedCoord::new(0.0, 10.0), 5, 10).unwrap(),
898        ];
899        let variogram = VariogramModel::new(0.01, 1.0, 30.0, VariogramType::Exponential).unwrap();
900        let prior = BinomialPrior::default();
901        let model = BinomialProjectedKrigingModel::new_with_prior(
902            obs.clone(),
903            variogram,
904            Anisotropy2D::isotropic(),
905            prior,
906        )
907        .unwrap()
908        .model;
909        let pred = model.predict(obs[0].coord()).unwrap();
910        let expected_logit = obs[0].smoothed_logit_with_prior(prior);
911        // Per-site logit observation noise: collocated prediction need not match local EB logit.
912        assert!(
913            (pred.logit_value - expected_logit).abs() < 0.35,
914            "logit err"
915        );
916        assert!(pred.prevalence > 0.0 && pred.prevalence < 1.0);
917        assert!(pred.variance >= 0.0);
918        assert!(pred.prevalence_variance >= 0.0);
919    }
920
921    #[test]
922    fn binomial_projected_anisotropy_pulls_predictions_along_major_axis() {
923        // Two stations along x-axis with low prevalence, two along y-axis with high.
924        let obs = vec![
925            ProjectedBinomialObservation::new(ProjectedCoord::new(-5.0, 0.0), 1, 10).unwrap(),
926            ProjectedBinomialObservation::new(ProjectedCoord::new(5.0, 0.0), 1, 10).unwrap(),
927            ProjectedBinomialObservation::new(ProjectedCoord::new(0.0, 5.0), 9, 10).unwrap(),
928            ProjectedBinomialObservation::new(ProjectedCoord::new(0.0, -5.0), 9, 10).unwrap(),
929        ];
930        let variogram = VariogramModel::new(0.01, 1.0, 20.0, VariogramType::Exponential).unwrap();
931        let iso =
932            BinomialProjectedKrigingModel::new(obs.clone(), variogram, Anisotropy2D::isotropic())
933                .unwrap()
934                .model;
935        let aniso = BinomialProjectedKrigingModel::new(
936            obs.clone(),
937            variogram,
938            Anisotropy2D::new(0.0, 0.2).unwrap(),
939        )
940        .unwrap()
941        .model;
942        let target = ProjectedCoord::new(0.0, 0.0);
943        let iso_pred = iso.predict(target).unwrap();
944        let aniso_pred = aniso.predict(target).unwrap();
945        // Anisotropy weights x-axis (low prevalence) more heavily.
946        assert!(aniso_pred.prevalence < iso_pred.prevalence);
947    }
948
949    #[test]
950    fn binomial_projected_rejects_too_few_observations() {
951        let one =
952            vec![ProjectedBinomialObservation::new(ProjectedCoord::new(0.0, 0.0), 1, 5).unwrap()];
953        let variogram = VariogramModel::new(0.01, 1.0, 10.0, VariogramType::Exponential).unwrap();
954        let r = BinomialProjectedKrigingModel::new(one, variogram, Anisotropy2D::isotropic());
955        assert!(r.is_err());
956    }
957}