Struct bayes_estimate::cholesky::UDU[][src]

pub struct UDU<N: RealField> {
    pub zero: N,
    pub one: N,
    pub minus_one: N,
}

Fields

zero: None: Nminus_one: N

Implementations

impl<N: RealField> UDU<N>[src]

pub fn new() -> UDU<N>[src]

pub fn UdUrcond<R: Dim, C: Dim>(UD: &MatrixMN<N, R, C>) -> N where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

Estimate the reciprocal condition number for inversion of the original PSD * matrix for which UD is the factor UdU’ or LdL’.

Additional columns are ignored. The rcond of the original matrix is simply the rcond of its d factor Using the d factor is fast and simple, and avoids computing any squares.

pub fn UCrcond<R: Dim, C: Dim>(&self, UC: &MatrixMN<N, R, C>) -> N where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

Estimate the reciprocal condition number for inversion of the original PSD matrix for which U is the factor UU’

The rcond of the original matrix is simply the square of the rcond of diagonal(UC).

pub fn UdUfactor_variant1<R: Dim, C: Dim>(
    &self,
    M: &mut MatrixMN<N, R, C>,
    n: usize
) -> N where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

In place Modified upper triangular Cholesky factor of a Positive definite or semi-definite matrix M.

Reference: A+G p.218 Upper Cholesky algorithm modified for UdU’

Numerical stability may not be as good as M(k,i) is updated from previous results. Algorithm has poor locality of reference and avoided for large matrices. Infinity values on the diagonal can be factorised.

Input: M, n=last column to be included in factorisation, Strict lower triangle of M is ignored in computation

Output: M as UdU’ factor

strict_upper_triangle(M) = strict_upper_triangle(U), /// diagonal(M) = d, strict_lower_triangle(M) is unmodified

Return: reciprocal condition number, -1 if negative, 0 if semi-definite (including zero)

pub fn UdUfactor_variant2<R: Dim, C: Dim>(
    &self,
    M: &mut MatrixMN<N, R, C>,
    n: usize
) -> N where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

In place modified upper triangular Cholesky factor of a Positive definite or semi-definite matrix M

Reference: A+G p.219 right side of table

Algorithm has good locality of reference and preferable for large matrices. Infinity values on the diagonal cannot be factorised

Input: M, n=last column to be included in factorisation, Strict lower triangle of M is ignored in computation

Output: M as UdU’ factor

strict_upper_triangle(M) = strict_upper_triangle(U), diagonal(M) = d, strict_lower_triangle(M) is unmodified

Return: reciprocal condition number, -1 if negative, 0 if semi-definite (including zero)

pub fn UCfactor_n<R: Dim, C: Dim>(
    &self,
    M: &mut MatrixMN<N, R, C>,
    n: usize
) -> N where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

In place upper triangular Cholesky factor of a Positive definite or semi-definite matrix M.

Reference: A+G p.218

Input: M, n=last std::size_t to be included in factorisation, Strict lower triangle of M is ignored in computation

Output: M as UC*UC’ factor, upper_triangle(M) = UC

Return: reciprocal condition number, -1 if negative, 0 if semi-definite (including zero)

pub fn UdUinverse<R: Dim, C: Dim>(&self, UD: &mut MatrixMN<N, R, C>) -> bool where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

In-place (destructive) inversion of diagonal and unit upper triangular matrices in UD.

BE VERY CAREFUL THIS IS NOT THE INVERSE OF UD.

Inversion on d and U is separate: inv(U)*inv(d)*inv(U’) = inv(U’dU) NOT EQUAL inv(UdU’)

Lower triangle of UD is ignored and unmodified, Only diagonal part d can be singular (zero elements), inverse is computed of all elements other then singular.

Reference: A+G p.223

Output: UD: inv(U), inv(d)

Return: singularity (of d), true iff d has a zero element

pub fn UTinverse<R: Dim, C: Dim>(&self, U: &mut MatrixMN<N, R, C>) -> bool where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

In-place (destructive) inversion of upper triangular matrix in U.

Output: U: inv(U)

Return: singularity (of U), true iff diagonal of U has a zero element

pub fn UdUrecompose_transpose<R: Dim, C: Dim>(M: &mut MatrixMN<N, R, C>) where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

In-place recomposition of Symmetric matrix from U’dU factor store in UD format.

Generally used for recomposing result of UdUinverse. Note definiteness of result depends purely on diagonal(M) i.e. if d is positive definite (>0) then result is positive definite

Reference: A+G p.223

In place computation uses simple structure of solution due to triangular zero elements. Defn: R = (U’ d) row i , C = U column j -> M(i,j) = R dot C, However M(i,j) only dependent R(k<=i), C(k<=j) due to zeros. Therefore in place multiple sequences such k < i <= j

Input: M - U’dU factorisation (UD format)

Output: M - U’dU recomposition (symmetric)

pub fn UdUrecompose<R: Dim, C: Dim>(M: &mut MatrixMN<N, R, C>) where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

In-place recomposition of Symmetric matrix from UdU’ factor store in UD format.

See UdUrecompose_transpose()

Input: M - UdU’ factorisation (UD format)

Output: M - UdU’ recomposition (symmetric)

pub fn Lzero<R: Dim, C: Dim>(&self, M: &mut MatrixMN<N, R, C>) where
    DefaultAllocator: Allocator<N, R, C>, 
[src]

Zero strict lower triangle of Matrix.

Auto Trait Implementations

impl<N> RefUnwindSafe for UDU<N> where
    N: RefUnwindSafe

impl<N> Send for UDU<N>

impl<N> Sync for UDU<N>

impl<N> Unpin for UDU<N> where
    N: Unpin

impl<N> UnwindSafe for UDU<N> where
    N: 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> Same<T> for T

type Output = T

Should always be Self

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

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>,