[][src]Struct petal_decomposition::Pca

pub struct Pca<A> where
    A: Scalar, 
{ /* fields omitted */ }

Principal component analysis.

This reduces the dimensionality of the input data using Singular Value Decomposition (SVD). The data is centered for each feature before applying SVD.

Examples

use petal_decomposition::Pca;

let x = ndarray::arr2(&[[0_f64, 0_f64], [1_f64, 1_f64], [2_f64, 2_f64]]);
let y = Pca::new(1).fit_transform(&x).unwrap();  // [-2_f64.sqrt(), 0_f64, 2_f64.sqrt()]
assert!((y[(0, 0)].abs() - 2_f64.sqrt()).abs() < 1e-8);
assert!(y[(1, 0)].abs() < 1e-8);
assert!((y[(2, 0)].abs() - 2_f64.sqrt()).abs() < 1e-8);

Implementations

impl<A> Pca<A> where
    A: Scalar + Lapack,
    A::Real: ScalarOperand
[src]

#[must_use]pub fn new(n_components: usize) -> Self[src]

Creates a PCA model with the given number of components.

pub fn components(&self) -> &Array2<A>[src]

Returns the principal axes in feature space.

pub fn mean(&self) -> &Array1<A>[src]

Returns the per-feature empirical mean.

pub fn n_components(&self) -> usize[src]

Returns the number of components.

pub fn singular_values(&self) -> &Array1<A::Real>[src]

Returns sigular values.

pub fn explained_variance_ratio(&self) -> Array1<A::Real>[src]

Returns the ratio of explained variance for each component.

pub fn fit<S>(
    &mut self,
    input: &ArrayBase<S, Ix2>
) -> Result<(), DecompositionError> where
    S: Data<Elem = A>, 
[src]

Fits the model with input.

Errors

  • DecompositionError::InvalidInput if any of the dimensions of input is less than the number of components.
  • DecompositionError::LinalgError if the underlying Singular Vector Decomposition routine fails.

pub fn transform<S>(
    &self,
    input: &ArrayBase<S, Ix2>
) -> Result<Array2<A>, DecompositionError> where
    S: Data<Elem = A>, 
[src]

Applies dimensionality reduction to input.

Errors

  • DecompositionError::InvalidInput if the number of features in input does not match that of the training data.

pub fn fit_transform<S>(
    &mut self,
    input: &ArrayBase<S, Ix2>
) -> Result<Array2<A>, DecompositionError> where
    S: Data<Elem = A>, 
[src]

Fits the model with input and apply the dimensionality reduction on input.

This is equivalent to calling both fit and transform for the same input, but more efficient.

Errors

Returns DecompositionError::LinalgError if the underlying Singular Vector Decomposition routine fails.

pub fn inverse_transform<S>(
    &self,
    input: &ArrayBase<S, Ix2>
) -> Result<Array2<A>, DecompositionError> where
    S: Data<Elem = A>, 
[src]

Transforms data back to its original space.

Errors

Returns DecompositionError::InvalidInput if the number of rows of input is different from that of the training data, or the number of columns of input is different from the number of components.

Auto Trait Implementations

impl<A> RefUnwindSafe for Pca<A> where
    A: RefUnwindSafe,
    <A as Scalar>::Real: RefUnwindSafe

impl<A> Send for Pca<A> where
    A: Send,
    <A as Scalar>::Real: Send

impl<A> Sync for Pca<A> where
    A: Sync,
    <A as Scalar>::Real: Sync

impl<A> Unpin for Pca<A> where
    <A as Scalar>::Real: Unpin

impl<A> UnwindSafe for Pca<A> where
    A: RefUnwindSafe,
    <A as Scalar>::Real: RefUnwindSafe + UnwindSafe

Blanket Implementations

impl<T> Any for T where
    T: 'static + ?Sized
[src]

impl<T> Borrow<T> for T where
    T: ?Sized
[src]

impl<T> BorrowMut<T> for T where
    T: ?Sized
[src]

impl<T> From<T> for T[src]

impl<T, U> Into<U> for T where
    U: From<T>, 
[src]

impl<T, U> TryFrom<U> for T where
    U: Into<T>, 
[src]

type Error = Infallible

The type returned in the event of a conversion error.

impl<T, U> TryInto<U> for T where
    U: TryFrom<T>, 
[src]

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

The type returned in the event of a conversion error.

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