CurveFit

Struct CurveFit 

Source
pub struct CurveFit<'data, B, T: Value = f64>
where B: Basis<T> + PolynomialDisplay<T>,
{ /* private fields */ }
Expand description

Represents a polynomial curve fit for a set of data points.

CurveFit computes a polynomial that best fits a given dataset using a specified polynomial basis (e.g., monomial, Chebyshev). It stores both the original data and the resulting coefficients.

§For beginners

Most users do not need to construct this directly. Use one of the specialized type aliases for common bases:

§How it works

  • Builds a basis matrix with shape [rows, k] where rows is the number of data points and k is the number of basis functions.
  • Forms a column vector b from the y values of the dataset.
  • Solves the linear system A * x = b using the SVD of the basis matrix. The solution x is the vector of polynomial coefficients.

§Type parameters

  • B: The basis type, implementing Basis<T>.
  • T: Numeric type (default f64) implementing Value.

§Example

let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = MonomialFit::new(data, 2).unwrap();
println!("Coefficients: {:?}", fit.coefficients());

Implementations§

Source§

impl<'data, T: Value, B> CurveFit<'data, B, T>
where B: Basis<T> + PolynomialDisplay<T>,

Source

pub fn to_owned(&self) -> CurveFit<'static, B, T>

Returns an owned version of this curve fit, with a full copy of the data.

Source

pub fn new(data: impl Into<Cow<'data, [(T, T)]>>, degree: usize) -> Result<Self>

Creates a new polynomial curve fit for the given data and degree.

You can also use CurveFit::new_auto to automatically select the best degree.

This constructor fits a polynomial to the provided (x, y) points using the chosen basis type B. For most users, B will be a Chebyshev or Monomial basis.

§Parameters
  • data: Slice of (x, y) points to fit.
  • degree: Desired polynomial degree.
§Returns

Returns Ok(Self) if the fit succeeds.

§Errors

Returns an Error in the following cases:

  • Error::NoData: data is empty.
  • Error::DegreeTooHigh: degree >= data.len().
  • Error::Algebra: the linear system could not be solved.
  • Error::CastFailed: a numeric value could not be cast to the target type.
§Behavior
  • Builds the basis matrix internally and fills each row using Basis::fill_matrix_row.
  • Computes x_range as the inclusive range of x values in the data.
  • Solves the linear system A * x = b to determine polynomial coefficients.
  • Uses SVD to solve the system for numerical stability.
    • Will limit iterations to between 1,000 and 10,000 based on problem size.
    • Prevents infinite loops by capping iterations.
    • And if you have not converged by 10,000 iterations, you probably won’t.
    • If you hit that, use an orthogonal basis like Chebyshev or Legendre for better stability.
§Warning

For datasets with >1,000,000 points, the basis matrix may be constructed in parallel. The normal equation is used to reduce the system size only if the data is stable (see below).

For a given problem of NxK, where N is the number of data points and K is the number of basis functions, if N > 1e6, and condition number of (X^T X) > 10^(digits - 8), then the normal equation will be used. where `digits` is the number of significant digits for type `T`.

We use Kahan summation to reduce numerical error when accumulating the partial results, as well as Tikhonov regularization to improve stability.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
println!("Coefficients: {:?}", fit.coefficients());
Source

pub fn new_weighted<F>( sample_fn: F, x_range: RangeInclusive<T>, degree: usize, ) -> Result<Self>
where B: OrthogonalBasis<T>, F: Fn(T) -> T,

Creates a new weighted polynomial curve fit using Gauss nodes.

This constructor fits a polynomial to the provided sampling function evaluated at Gauss nodes, using weights associated with those nodes.

§Parameters
  • sample_fn: Function that takes an x value and returns the corresponding y value.
  • x_range: Inclusive range of x values to consider for fitting.
  • degree: Desired polynomial degree.
§Returns

Returns Ok(Self) if the fit succeeds.

§Errors

Returns an Error in the following cases:

  • Error::DegreeTooHigh: degree is invalid for the chosen basis.
  • Error::Algebra: the linear system could not be solved.
  • Error::CastFailed: a numeric value could not be cast to the target type.
Source

pub fn new_auto( data: impl Into<Cow<'data, [(T, T)]>>, max_degree: impl Into<DegreeBound>, method: &impl ModelScoreProvider, ) -> Result<Self>

Automatically selects the best polynomial degree and creates a curve fit.

This function fits polynomials of increasing degree to the provided dataset and selects the “best” degree according to the specified scoring method.

§Parameters
  • data: Slice of (x, y) points to fit.
  • method: crate::score to evaluate model quality.
    • AIC: Akaike Information Criterion (uses AICc if n/k < 4)
    • BIC: Bayesian Information Criterion
§Choosing a scoring method
  • Consider the size of your dataset: If you have a small dataset, prefer AIC as it penalizes complexity more gently.
  • If your dataset is large, BIC may be more appropriate as it imposes a harsher penalty on complexity.
§Returns

Returns Ok(Self) with the fit at the optimal degree.

§Errors

Returns Error if:

  • data is empty (Error::NoData)
  • A numeric value cannot be represented in the target type (Error::CastFailed)
§Behavior
  • Starts with degree 0 and iteratively fits higher degrees up to data.len() - 1.
  • Evaluates each fit using model_score(method).
  • Selects best model, where the score is within 2 of the minimum score
    • As per Burnham and Anderson, models within Δ ≤ 2 are statistically indistinguishable from the best model.
  • Uses SVD to solve the system for numerical stability.
    • Will limit iterations to between 1,000 and 10,000 based on problem size.
    • Prevents infinite loops by capping iterations.
    • And if you have not converged by 10,000 iterations, you probably won’t.
    • If you hit that, use an orthogonal basis like Chebyshev or Legendre for better stability.
§Warning

For datasets with >1,000,000 points, the basis matrix may be constructed in parallel. The normal equation is used to reduce the system size only if the data is stable (see below).

For a given problem of NxK, where N is the number of data points and K is the number of basis functions, if N > 1e6, and condition number of (X^T X) > 10^(digits - 8), then the normal equation will be used.

We use Kahan summation to reduce numerical error when accumulating the partial results, as well as Tikhonov regularization to improve stability.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new_auto(data, DegreeBound::Relaxed, &Aic).unwrap();
println!("Selected degree: {}", fit.degree());
Source

pub fn new_kfold_cross_validated( data: impl Into<Cow<'data, [(T, T)]>>, strategy: CvStrategy, max_degree: impl Into<DegreeBound>, method: &impl ModelScoreProvider, ) -> Result<Self>

Creates a new polynomial curve fit using K-fold cross-validation to select the best degree.

This function splits the dataset into folds subsets, using each subset as a validation set while training on the remaining data. It evaluates polynomial fits of increasing degree and selects the best degree based on the specified scoring method.

This method helps prevent overfitting by ensuring that the selected model generalizes well to unseen data, and is particularly useful for small datasets or those with very extreme outliers.

§Parameters
  • data: Slice of (x, y) points to fit.
  • strategy: Cross-validation strategy defining the number of folds.
  • max_degree: Maximum polynomial degree to consider.
  • method: ModelScoreProvider to evaluate model quality.
§Choosing a scoring method
  • Consider the size of your dataset: If you have a small dataset, prefer AIC as it penalizes complexity more gently.
  • If your dataset is large, BIC may be more appropriate as it imposes a harsher penalty on complexity.
§Choosing number of folds (strategy)

The number of folds (k) determines how many subsets the data is split into for cross-validation. The choice of k can impact the bias-variance trade-off:

  • Bias is how much your model’s predictions differ from the actual values on average

  • Variance is how much your model’s predictions change when you use different training data

  • Choose [CvStrategy::MinimizeBias] (e.g., k=5) with larger datasets to reduce bias. Prevents a model from being too simple to capture data patterns

  • Choose [CvStrategy::MinimizeVariance] (e.g., k=10) with smaller datasets to reduce variance. Helps ensure the model generalizes well to unseen data.

  • Choose [CvStrategy::Balanced] (e.g., k=7) for a compromise between bias and variance.

  • Choose [CvStrategy::LeaveOneOut] (k=n) for very small datasets, but be aware of the high computational cost.

  • Choose [CvStrategy::Custom(k)] to specify a custom number of folds, ensuring k >= 2 and k <= n.

§Returns

Returns Ok(Self) with the fit at the optimal degree.

§Errors

Returns Error if:

  • data is empty or folds < 2 (Error::NoData)
  • A numeric value cannot be represented in the target type (Error::CastFailed)
  • No valid model could be fitted (Error::NoModel)
§Warning

For datasets with >1,000,000 points, the basis matrix may be constructed in parallel. The normal equation is used to reduce the system size only if the data is stable (see below).

For a given problem of NxK, where N is the number of data points and K is the number of basis functions, if N > 1e6, and condition number of (X^T X) > 10^(digits - 8), then the normal equation will be used.

We use Kahan summation to reduce numerical error when accumulating the partial results, as well as Tikhonov regularization to improve stability.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new_kfold_cross_validated(data, CvStrategy::LeaveOneOut, DegreeBound::Relaxed, &Aic).unwrap();
println!("Selected degree: {}", fit.degree());
Source

pub fn prune_insignificant( &mut self, confidence: Confidence, ) -> Result<Vec<(usize, T)>>

Prunes coefficients that are statistically insignificant based on a t-test.

§Parameters
  • confidence: Confidence level for determining significance (e.g., P95, P99)
§Errors

Returns an error if the covariance matrix cannot be computed.

§Returns

A vector of (index, coefficient) for all pruned coefficients.

§Notes
  • Modifies self in-place, zeroing out insignificant coefficients.
  • Uses the standard errors derived from the covariance matrix.
  • Ignores coefficients whose absolute value is smaller than T::epsilon().
Source

pub fn covariance(&self) -> Result<CurveFitCovariance<'_, '_, B, T>>

Computes the covariance matrix and related statistics for this curve fit.

The returned CurveFitCovariance provides:

  • Covariance matrix of the fitted coefficients.
  • Standard errors of the coefficients.
  • Prediction variance at a specific x.
  • Confidence intervals for predictions.
§Returns
  • Ok(CurveFitCovariance<'_, B, T>) on success.
§Errors

Returns Err(Error::Algebra) if (Xᵀ X) is singular or nearly singular, i.e., the pseudo-inverse cannot be computed. Causes include too few data points relative to parameters or collinear/linearly dependent basis functions.

§Example
let cov = model.covariance().unwrap();
let se = cov.coefficient_standard_errors();
let band = cov.confidence_band(1.0, Confidence::P95, None).unwrap();
println!("Predicted CI at x=1: {} - {}", band.min(), band.max());
Source

pub fn critical_points(&self) -> Result<Vec<CriticalPoint<T>>>

Finds the critical points (where the derivative is zero) of a polynomial in this basis.

This corresponds to the polynomial’s local minima and maxima (The x values where curvature changes).

Technical Details

The critical points are found by solving the equation f'(x) = 0, where f'(x) is the derivative of the polynomial.

This is done with by finding the eigenvalues of the companion matrix of the derivative polynomial.

§Returns

A vector of x values where the critical points occur.

§Requirements
§Errors

Returns an error if the critical points cannot be found.

§Example
let critical_points = model.critical_points().unwrap();
Source

pub fn area_under_curve( &self, x_min: T, x_max: T, constant: Option<T>, ) -> Result<T>
where B: IntegralBasis<T>,

Computes the definite integral (area under the curve) of the fitted polynomial between x_min and x_max.

Technical Details

The area under the curve is computed using the definite integral of the polynomial between the specified bounds:

Area = ∫[x_min to x_max] f(x) dx = F(x_max) - F(x_min)
§Parameters
  • x_min: Lower bound of integration.
  • x_max: Upper bound of integration.
  • constant: Constant of integration (value at x = 0) for the indefinite integral.
§Requirements
§Returns
  • Ok(T): The computed area under the curve between x_min and x_max.
  • Err: If computing the integral fails (e.g., basis cannot compute integral coefficients).
§Errors

If the basis cannot compute the integral coefficients, an error is returned.

§Example
let area = model.area_under_curve(0.0, 3.0, None).unwrap();
println!("Area under curve: {}", area);
Source

pub fn monotonicity_violations(&self) -> Result<Vec<T>>
where B: DifferentialBasis<T>,

Returns the X-values where the function is not monotone (i.e., where the derivative changes sign).

§Errors

Returns an error if the derivative cannot be computed.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = MonomialFit::new(data, 2).unwrap();
let violations = fit.monotonicity_violations().unwrap();
Source

pub fn model_score(&self, method: &impl ModelScoreProvider) -> T

Computes the quality score of the polynomial fit using the specified method.

This evaluates how well the fitted polynomial represents the data, taking into account both the fit error and model complexity.

§Parameters
  • method: ModelScoreProvider to use for scoring.
    • AIC: Akaike Information Criterion (uses AICc if n/k < 4)
    • BIC: Bayesian Information Criterion
§Returns

The score as a numeric value (T). Lower scores indicate better models.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let score = fit.model_score(&Aic);
println!("Model score: {}", score);
Source

pub fn residuals(&self) -> Vec<(T, T)>

Computes the residuals of the fit.

Residuals are the differences between the observed y values and the predicted y values from the fitted polynomial. They provide insight into the fit quality and can be used for diagnostic purposes.

Technical Details

residual_i = y_i - f(x_i)
where
  y_i = observed value, f(x_i) = predicted value from the polynomial at x_i
§Returns

A vector of residuals, where each element corresponds to a data point.

Source

pub fn filtered_residuals(&self) -> Vec<(T, T)>

Computes the residuals of the fit, filtering out small residuals likely due to floating point noise.

This form can help minimize the impact of floating point precision and rounding.

Technical Details

max_val = max(|y_i|, |f(x_i)|, 1)
epsilon = machine_epsilon * sqrt(n) * max_val
residual_i = y_i - f(x_i) if |y_i - f(x_i)| >= epsilon else 0
where
  y_i = observed value, f(x_i) = predicted value from the polynomial at x_i, n = number of data points
§Returns

A vector of scaled residuals, where each element corresponds to a data point.

Source

pub fn residual_variance(&self) -> T

Computes the residual variance of the model’s predictions.

See statistics::residual_variance.

Residual variance is the unbiased estimate of the variance of the errors (σ²) after fitting a model. It’s used for confidence intervals and covariance estimates of the fitted parameters.

Source

pub fn mean_squared_error(&self) -> T

Computes the mean squared error (MSE) of this fit against its source data.

See statistics::mean_squared_error.

MSE measures the average squared difference between the observed and predicted values. Lower values indicate a better fit.

Source

pub fn root_mean_squared_error(&self) -> T

Computes the root mean squared error (RMSE) of this fit against its source data.

See statistics::root_mean_squared_error.

RMSE is the square root of the MSE, giving error in the same units as the original data. Lower values indicate a closer fit.

Source

pub fn mean_absolute_error(&self) -> T

Computes the mean absolute error (MAE) of this fit against its source data.

See statistics::mean_absolute_error.

MAE measures the average absolute difference between observed and predicted values. Lower values indicate a better fit.

Source

pub fn r_squared(&self, data: Option<&[(T, T)]>) -> T

Calculates the R-squared value for the model compared to provided data.

If not provided, uses the original data used to create the fit.

That way you can test the fit against a different dataset if desired.

R-squared is a statistical measure of how well the polynomial explains the variance in the data. Values closer to 1 indicate a better fit.

§Parameters
  • data: A slice of (x, y) pairs to compare against the polynomial fit.

See statistics::r_squared for more details.

§Returns

The R-squared value as type T.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let r2 = fit.r_squared(None);
println!("R² = {}", r2);
Source

pub fn robust_r_squared(&self, data: Option<&[(T, T)]>) -> T

Uses huber loss to compute a robust R-squared value.

  • More robust to outliers than traditional R².
  • Values closer to 1 indicate a better fit.

If Data not provided, uses the original data used to create the fit.

That way you can test the fit against a different dataset if desired.

Technical Details

R²_robust = 1 - (Σ huber_loss(y_i - y_fit_i, delta)) / (Σ (y_i - y_fit_i)²)
where
  huber_loss(r, delta) = { 0.5 * r²                    if |r| ≤ delta
                         { delta * (|r| - 0.5 * delta) if |r| > delta
 delta = 1.345 * MAD
 MAD = median( |y_i - y_fit_i| )
 where
   y_i = observed values, y_fit_i = predicted values
§Type Parameters
  • T: A numeric type implementing the Value trait.
§Parameters
  • data: A slice of (x, y) pairs to compare against the polynomial fit.
§Returns

The robust R² as a T.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let r2 = fit.robust_r_squared(None);
println!("R² = {}", r2);
Source

pub fn adjusted_r_squared(&self, data: Option<&[(T, T)]>) -> T

Computes the adjusted R-squared value.

Adjusted R² accounts for the number of predictors in a model, penalizing overly complex models. Use it to compare models of different degrees.

If Data not provided, uses the original data used to create the fit.

That way you can test the fit against a different dataset if desired.

Technical Details

R²_adj = R² - (1 - R²) * (k / (n - k))
where
  n = number of observations, k = number of model parameters

[r_squared] is used to compute R²

§Type Parameters
  • T: A numeric type implementing the Value trait.
§Parameters
  • data: A slice of (x, y) pairs to compare against
§Returns

The adjusted R² as a T.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let r2 = fit.adjusted_r_squared(None);
println!("R² = {}", r2);
Source

pub fn r_squared_against<C>(&self, function: &Polynomial<'_, C, T>) -> T
where C: Basis<T> + PolynomialDisplay<T>,

Calculates the R-squared value for the model compared to provided function.

R-squared is a statistical measure of how well the polynomial explains the variance in the data. Values closer to 1 indicate a better fit.

§Parameters
  • data: A slice of (x, y) pairs to compare against the polynomial fit.

See statistics::r_squared for more details.

§Returns

The R-squared value as type T.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let target = MonomialPolynomial::borrowed(&[1.0, 2.0, 1.0]);
let r2 = fit.r_squared_against(&target);
println!("R² vs target polynomial = {}", r2);
Source

pub fn folded_rmse(&self, strategy: CvStrategy) -> UncertainValue<T>

Computes the folded root mean squared error (RMSE) with uncertainty estimation.

This is a measure of how far the predicted values are from the actual values, on average.

Doing it over K-fold cross-validation (Splitting the data into K subsets, and leaving one out for testing each time) tells us how the model performance changes as the data varies.

When to use each strategy:

  • MinimizeBias: When the dataset is small and you want to avoid underfitting. Prevents a model from being too simple to capture data patterns.
  • MinimizeVariance: When the dataset is large and you want to avoid overfitting. Helps ensure the model generalizes well to unseen data.
  • Balanced: When you want a good trade-off between bias and variance, suitable for moderately sized datasets or when unsure.
  • Custom: Specify your own number of folds (k) based on domain knowledge or specific requirements. Use with caution

The returned statistics::UncertainValue<T> includes both the mean RMSE and its uncertainty (standard deviation across folds).

You can use statistics::UncertainValue::confidence_band to get confidence intervals for the RMSE.

§Parameters
  • strategy: The cross-validation strategy to use (e.g., K-Fold with K=5).
§Returns

An statistics::UncertainValue<T> representing the folded RMSE with uncertainty.

Source

pub fn degree(&self) -> usize

Returns the degree of the polynomial.

The number of actual components, or basis functions, in the expression of a degree is defined by the basis.

That number is called k. For most basis choices, k = degree + 1.

Source

pub fn coefficients(&self) -> &[T]

Returns a reference to the polynomial’s coefficients.

The index of each coefficient the jth basis function.

For example in a monomial expression y(x) = 2x^2 - 3x + 1; coefficients = [1.0, -3.0, 2.0]

Technical Details

Formally, for each coefficient j, and the jth basis function B_j(x), the relationship is:

y(x) = Σ (c_j * B_j(x))
Source

pub fn data(&self) -> &[(T, T)]

Returns a reference to the data points used for fitting.

Each element is a (x, y) tuple representing a data point.

Source

pub fn x_range(&self) -> RangeInclusive<T>

Returns the inclusive range of x-values in the dataset.

Source

pub fn y_range(&self) -> RangeInclusive<T>

Returns the inclusive range of y-values in the dataset.

This is computed dynamically from the stored data points. Use sparingly

Source

pub fn y(&self, x: T) -> Result<T>

Evaluates the polynomial at a given x-value.

Technical Details

Given Basis::k coefficients and basis functions, and for each pair of coefficients c_j and basis function B_j(x), this function returns:

y(x) = Σ (c_j * B_j(x))
§Parameters
  • x: The point at which to evaluate the polynomial.
§Returns

The corresponding y-value as T if x is within the valid range.

§Errors

Returns Error::DataRange if x is outside the original data bounds.

§Notes
  • Polynomial fits are generally only stable within the x-range used for fitting.
  • To evaluate outside the original bounds, use CurveFit::as_polynomial to get a pure polynomial function that ignores the original x-range.
§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let y = fit.y(1.0).unwrap();
println!("y(1.0) = {}", y);
Source

pub fn solution(&self) -> Vec<(T, T)>

Returns the fitted y-values corresponding to the original x-values.

This produces a vector of (x, y) pairs for the same x-values used in the source data. It is guaranteed to succeed because all x-values are within the curve’s valid range.

§Notes
  • Useful for quickly plotting or analyzing the fitted curve against the original data points.
  • The method internally calls CurveFit::y but is infallible because it only evaluates x-values within the valid range.
§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let points = fit.solution();
for (x, y) in points {
    println!("x = {}, y = {}", x, y);
}
Source

pub fn solve(&self, x: impl IntoIterator<Item = T>) -> Result<Vec<(T, T)>>

Evaluates the curve at multiple x-values.

Technical Details

Given Basis::k coefficients and basis functions, and for each pair of coefficients c_j and basis function B_j(x), this function returns:

y(x) = Σ (c_j * B_j(x))
§Parameters
  • x: An iterator of x-values to evaluate.
§Returns

A vector of (x, y) pairs corresponding to each input x-value.

§Errors

Returns Error::DataRange if any x-value is outside the original data range.

§Notes
  • Curve fits are generally only stable within the x-range used for fitting.
  • To evaluate outside the original bounds, use CurveFit::as_polynomial to get a pure polynomial function that ignores the original x-range.
§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let points = fit.solve([0.0, 0.5, 1.0]).unwrap();
for (x, y) in points {
    println!("x = {}, y = {}", x, y);
}
Source

pub fn solve_range( &self, range: RangeInclusive<T>, step: T, ) -> Result<Vec<(T, T)>>

Evaluates the curve at evenly spaced points over a range.

Technical Details

Given Basis::k coefficients and basis functions, and for each pair of coefficients c_j and basis function B_j(x), this function returns:

y(x) = Σ (c_j * B_j(x))
§Parameters
  • range: The start and end x-values to evaluate.
  • step: The increment between points.
§Returns

A vector of (x, y) pairs corresponding to each x-value in the range.

§Errors

Returns Error::DataRange if any x-value is outside the original data range.

§Notes
  • Curve fits are only stable within the x-range used for fitting.
  • To evaluate outside the original bounds, use CurveFit::as_polynomial to get a pure polynomial function.
§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let points = fit.solve_range(0.0..=2.0, 0.5).unwrap();
for (x, y) in points {
    println!("x = {}, y = {}", x, y);
}
Source

pub fn as_polynomial(&self) -> &Polynomial<'_, B, T>

Returns a pure polynomial representation of the curve fit.

This allows evaluation of the polynomial at any x-value, without restriction to the original data range. Unlike CurveFit::y or CurveFit::solve, this does not perform range checks, so use with caution outside the fit’s stable region.

The Polynomial form is considered a canonical function, not a fit estimate.

§Returns

A reference to the Polynomial that this fit uses internally.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let poly = fit.as_polynomial();
let y = poly.y(10.0); // can evaluate outside original x-range
Source

pub fn into_polynomial(self) -> Polynomial<'static, B, T>

Returns a pure polynomial representation of the curve fit.

This is primarily a memory optimization in practice; You become responsible for maintaining the stability around the x-bounds, but the copy of the original data, and the Vandermonde matrix are dropped.

This allows evaluation of the polynomial at any x-value, without restriction to the original data range. Unlike CurveFit::y or CurveFit::solve, this does not perform range checks, so use with caution outside the fit’s stable region.

The Polynomial form is considered a canonical function, not a fit estimate.

§Returns

The Polynomial that this fit uses internally

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let poly = fit.as_polynomial();
let y = poly.y(10.0); // can evaluate outside original x-range
Source

pub fn as_monomial(&self) -> Result<MonomialPolynomial<'static, T>>
where B: IntoMonomialBasis<T>,

Converts the curve fit into a monomial polynomial.

This produces a MonomialPolynomial representation of the curve, which uses the standard monomial basis 1, x, x^2, ….

§Returns

A monomial polynomial with owned coefficients.

§Errors

Returns an error if the current basis cannot be converted to monomial form. This requires that the basis implements IntoMonomialBasis.

§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
let mono_poly = fit.as_monomial().unwrap();
let y = mono_poly.y(1.5);
Source

pub fn coefficient_energies(&self) -> Result<Vec<T>>
where B: OrthogonalBasis<T>,

Computes the energy contribution of each coefficient in an orthogonal basis.

This is a measure of how much each basis function contributes to the resulting polynomial.

It can be useful for understanding the significance of each term

Technical Details

For an orthogonal basis, the energy contribution of each coefficient is calculated as:

E_j = c_j^2 * N_j

where:

  • E_j is the energy contribution of the jth coefficient.
  • c_j is the jth coefficient.
  • N_j is the normalization factor for the jth basis function, provided by the basis.
§Returns

A vector of energy contributions for each coefficient.

§Errors

Returns an error if the series is not orthogonal, like with integrated Fourier series.

Source

pub fn smoothness(&self) -> Result<T>
where B: OrthogonalBasis<T>,

Computes a smoothness metric for the polynomial.

This metric quantifies how “smooth” the polynomial is, with lower values indicating smoother curves.

Technical Details

The smoothness is calculated as a weighted average of the coefficient energies, where higher-degree coefficients are penalized more heavily. The formula used is:

Smoothness = (Σ (k^2 * E_k)) / (Σ E_k)

where:

  • k is the degree of the basis function.
  • E_k is the energy contribution of the k-th coefficient.
§Returns

A smoothness value, where lower values indicate a smoother polynomial.

§Errors

Returns an error if the series is not orthogonal, like with integrated Fourier series.

Source

pub fn spectral_energy_filter(&mut self) -> Result<()>
where B: OrthogonalBasis<T>,

Applies a spectral energy filter to the series.

This uses the properties of a orthogonal basis to de-noise the polynomial by removing higher-degree terms that contribute little to the overall energy. Terms are split into “signal” and “noise” based on their energy contributions, and the polynomial is truncated to only include the signal components.

Remaining terms are smoothly attenuated to prevent ringing artifacts from a hard cutoff.

Technical Details

The energy of each coefficient is calculated using the formula:

E_j = c_j^2 * N_j

where:

  • E_j is the energy contribution of the jth coefficient.
  • c_j is the jth coefficient.
  • N_j is the normalization factor for the jth basis function, provided by the basis.

Generalized Cross-Validation (GCV) is used to determine the optimal cutoff degree K that minimizes the prediction error using: GCV(K) = (suffix[0] - suffix[K]) / K^2, where suffix is the suffix sum of energies.

A Lanczos Sigma filter with p=1 is applied to smoothly attenuate coefficients up to the cutoff degree, reducing Gibbs ringing artifacts.

§Notes
  • This method modifies the polynomial in place.
§Errors

Returns an error if the basis is not orthogonal. This can be checked with Polynomial::is_orthogonal.

Source

pub fn equation(&self) -> String

Returns a human-readable string of the polynomial equation.

The output shows the polynomial in standard mathematical notation, for example:

y = 1.0x^3 + 2.0x^2 + 3.0x + 4.0
§Notes
  • Requires the basis to implement PolynomialDisplay for formatting.
  • This operation is infallible and guaranteed to succeed, hence no error return.
§Example
let data = &[(0.0, 1.0), (1.0, 3.0), (2.0, 7.0)];
let fit = ChebyshevFit::new(data, 2).unwrap();
println!("{}", fit.equation());
Source

pub fn properties(&self) -> FitProperties<T>

Returns the properties of the curve fit.

This is a comprehensive summary of the fit’s characteristics.

Trait Implementations§

Source§

impl<B, T> AsPlottingElement<T> for CurveFit<'_, B, T>
where B: Basis<T> + PolynomialDisplay<T>, T: Value,

Available on crate feature plotting only.
Source§

fn as_plotting_element( &self, _: &[T], confidence: Confidence, noise_tolerance: Option<Tolerance<T>>, ) -> PlottingElement<T>

Converts this to a plotting element
Source§

impl<B, T: Value> AsRef<Polynomial<'_, B, T>> for CurveFit<'_, B, T>
where B: Basis<T> + PolynomialDisplay<T>,

Source§

fn as_ref(&self) -> &Polynomial<'static, B, T>

Converts this type into a shared reference of the (usually inferred) input type.
Source§

impl<'data, B, T: Clone + Value> Clone for CurveFit<'data, B, T>
where B: Basis<T> + PolynomialDisplay<T> + Clone,

Source§

fn clone(&self) -> CurveFit<'data, B, T>

Returns a duplicate of the value. Read more
1.0.0 · Source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
Source§

impl<'data, B, T: Debug + Value> Debug for CurveFit<'data, B, T>
where B: Basis<T> + PolynomialDisplay<T> + Debug,

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

impl<T: Value, B> Display for CurveFit<'_, B, T>
where B: Basis<T> + PolynomialDisplay<T>,

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
Source§

impl<'data, B, T: PartialEq + Value> PartialEq for CurveFit<'data, B, T>
where B: Basis<T> + PolynomialDisplay<T> + PartialEq,

Source§

fn eq(&self, other: &CurveFit<'data, B, T>) -> bool

Tests for self and other values to be equal, and is used by ==.
1.0.0 · Source§

fn ne(&self, other: &Rhs) -> bool

Tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
Source§

impl<'data, B, T: Value> StructuralPartialEq for CurveFit<'data, B, T>
where B: Basis<T> + PolynomialDisplay<T>,

Auto Trait Implementations§

§

impl<'data, B, T> Freeze for CurveFit<'data, B, T>
where T: Freeze, B: Freeze,

§

impl<'data, B, T> RefUnwindSafe for CurveFit<'data, B, T>

§

impl<'data, B, T> Send for CurveFit<'data, B, T>

§

impl<'data, B, T> Sync for CurveFit<'data, B, T>

§

impl<'data, B, T> Unpin for CurveFit<'data, B, T>
where T: Unpin, B: Unpin,

§

impl<'data, B, T> UnwindSafe for CurveFit<'data, B, T>

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dest: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dest. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> IntoEither for T

Source§

fn into_either(self, into_left: bool) -> Either<Self, Self>

Converts self into a Left variant of Either<Self, Self> if into_left is true. Converts self into a Right variant of Either<Self, Self> otherwise. Read more
Source§

fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
where F: FnOnce(&Self) -> bool,

Converts self into a Left variant of Either<Self, Self> if into_left(&self) returns true. Converts self into a Right variant of Either<Self, Self> otherwise. Read more
Source§

impl<T> Pointable for T

Source§

const ALIGN: usize

The alignment of pointer.
Source§

type Init = T

The type for initializers.
Source§

unsafe fn init(init: <T as Pointable>::Init) -> usize

Initializes a with the given initializer. Read more
Source§

unsafe fn deref<'a>(ptr: usize) -> &'a T

Dereferences the given pointer. Read more
Source§

unsafe fn deref_mut<'a>(ptr: usize) -> &'a mut T

Mutably dereferences the given pointer. Read more
Source§

unsafe fn drop(ptr: usize)

Drops the object pointed to by the given pointer. Read more
Source§

impl<T> Same for T

Source§

type Output = T

Should always be Self
Source§

impl<SS, SP> SupersetOf<SS> for SP
where SS: SubsetOf<SP>,

Source§

fn to_subset(&self) -> Option<SS>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
Source§

fn is_in_subset(&self) -> bool

Checks if self is actually part of its subset T (and can be converted to it).
Source§

fn to_subset_unchecked(&self) -> SS

Use with care! Same as self.to_subset but without any property checks. Always succeeds.
Source§

fn from_subset(element: &SS) -> SP

The inclusion map: converts self to the equivalent element of its superset.
Source§

impl<T> ToOwned for T
where T: Clone,

Source§

type Owned = T

The resulting type after obtaining ownership.
Source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
Source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
Source§

impl<T> ToString for T
where T: Display + ?Sized,

Source§

fn to_string(&self) -> String

Converts the given value to a String. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
Source§

impl<V, T> VZip<V> for T
where V: MultiLane<T>,

Source§

fn vzip(self) -> V

Source§

impl<E, T> WithTypeFrom<T> for E
where E: AsPlottingElement<T>, T: Value,

Source§

fn options_with_type_from(&self) -> PlotOptions<T>

Available on crate feature plotting only.
Get default plot options for this type
Source§

fn plot_with_type_from<P>( &self, root: &<P as PlotBackend>::Root, options: PlotOptions<T>, ) -> Result<Plot<P, T>, <P as PlotBackend>::Error>
where P: PlotBackend,

Available on crate feature plotting only.
Create a new plot with this as the primary function Read more
Source§

fn add_to_plot_with_type_from<P>( &self, plot: &mut Plot<P, T>, ) -> Result<(), <P as PlotBackend>::Error>
where P: PlotBackend,

Available on crate feature plotting only.
Adds this element to the given plot Read more
Source§

impl<T> Scalar for T
where T: 'static + Clone + PartialEq + Debug,