ReFRACtor
|
Functions | |
blitz::Array< double, 2 > | FullPhysics::cholesky_decomposition (const blitz::Array< double, 2 > &A) |
This calculates the Cholesky Decompostion of A so that A = L L^T. More... | |
template<class T > | |
void | FullPhysics::fortran_resize (blitz::Array< T, 1 > &M, int sz1) |
If matrix is in Fortran order, resize it. More... | |
template<class T > | |
void | FullPhysics::fortran_resize (blitz::Array< T, 2 > &M, int sz1, int sz2) |
If matrix is in Fortran order, resize it. More... | |
template<class T > | |
void | FullPhysics::fortran_resize (blitz::Array< T, 3 > &M, int sz1, int sz2, int sz3) |
If matrix is in Fortran order, resize it. More... | |
blitz::Array< double, 2 > | FullPhysics::generalized_inverse (const blitz::Array< double, 2 > &A, double Rcond=std::numeric_limits< double >::epsilon()) |
This returns the generalized inverse of the given matrix. More... | |
blitz::Array< double, 2 > | FullPhysics::generalized_inverse (const blitz::Array< double, 2 > &A, const blitz::Array< bool, 1 > &Zero_unused, double Rcond=std::numeric_limits< double >::epsilon()) |
This returns the generalized inverse of the given matrix. More... | |
blitz::Array< double, 2 > | FullPhysics::multifit_covar (const blitz::Array< double, 2 > &J, double Eps_rel=std::numeric_limits< double >::epsilon()) |
bool | FullPhysics::multifit_test_delta (const blitz::Array< double, 1 > &dx, const blitz::Array< double, 1 > &x, double Eps_abs=std::numeric_limits< double >::epsilon(), double Eps_rel=std::numeric_limits< double >::epsilon()) |
bool | FullPhysics::multifit_test_gradient (const blitz::Array< double, 1 > &g, double Eps_abs=std::numeric_limits< double >::epsilon()) |
blitz::Array< double, 1 > | FullPhysics::solve_least_squares (const blitz::Array< double, 2 > &A, const blitz::Array< double, 1 > &B, double Rcond=1e-12) |
This solves the least squares system A*x = b, returning the minimum norm solution. More... | |
blitz::Array< double, 1 > | FullPhysics::solve_least_squares_qr (const blitz::Array< double, 2 > &A, const blitz::Array< double, 1 > &B) |
This finds the least squares solution to the overdetermined system A x = b where the matrix A has more rows than columns. More... | |
double | FullPhysics::sqr (double x) |
void | FullPhysics::svd (const blitz::Array< double, 2 > &A, blitz::Array< double, 1 > &S, blitz::Array< double, 2 > &U, blitz::Array< double, 2 > &VT) |
Compute the SVD decomposition of a matrix A = U * SIGMA * V^T. More... | |
template<class T , int D> | |
blitz::Array< T, D > | FullPhysics::to_c_order (const blitz::Array< T, D > &In) |
Ensure that a given blitz::Array is contiguous, not reversed, and in C RowMajorArray format. More... | |
template<class T , int D> | |
blitz::Array< T, D > | FullPhysics::to_c_order_const (const blitz::Array< T, D > &In) |
Ensure that a given blitz::Array is contiguous, not reversed, and in C RowMajorArray format. More... | |
template<class T , int D> | |
blitz::Array< T, D > | FullPhysics::to_fortran (const blitz::Array< T, D > &In) |
Ensure that a given blitz::Array is contiguous, not reversed, and in fortran ColumnMajorArray format. More... | |
template<class T , int D> | |
blitz::Array< T, D > | FullPhysics::to_fortran_const (const blitz::Array< T, D > &In) |
Ensure that a given blitz::Array is contiguous, not reversed, and in fortran ColumnMajorArray format. More... | |
blitz::Array< double, 2 > FullPhysics::cholesky_decomposition | ( | const blitz::Array< double, 2 > & | A | ) |
This calculates the Cholesky Decompostion of A so that A = L L^T.
This returns L. A must be symmetric positive definite.
Definition at line 212 of file linear_algebra.cc.
void FullPhysics::fortran_resize | ( | blitz::Array< T, 1 > & | M, |
int | sz1 | ||
) |
If matrix is in Fortran order, resize it.
Otherwise, set it to a Fortran matrix of the given size.
Definition at line 91 of file linear_algebra.h.
void FullPhysics::fortran_resize | ( | blitz::Array< T, 2 > & | M, |
int | sz1, | ||
int | sz2 | ||
) |
If matrix is in Fortran order, resize it.
Otherwise, set it to a Fortran matrix of the given size.
Definition at line 110 of file linear_algebra.h.
void FullPhysics::fortran_resize | ( | blitz::Array< T, 3 > & | M, |
int | sz1, | ||
int | sz2, | ||
int | sz3 | ||
) |
If matrix is in Fortran order, resize it.
Otherwise, set it to a Fortran matrix of the given size.
Definition at line 129 of file linear_algebra.h.
Array< double, 2 > FullPhysics::generalized_inverse | ( | const blitz::Array< double, 2 > & | A, |
double | Rcond = std::numeric_limits<double>::epsilon() |
||
) |
This returns the generalized inverse of the given matrix.
This is the same as the inverse, except that if A is rank deficient we set the SVD value to 0.
A | Matrix to get inverse of |
Rcond | We set singular values < Rcond * max singular value to 0. |
Definition at line 136 of file linear_algebra.cc.
Array< double, 2 > FullPhysics::generalized_inverse | ( | const blitz::Array< double, 2 > & | A, |
const blitz::Array< bool, 1 > & | Zero_unused, | ||
double | Rcond = std::numeric_limits<double>::epsilon() |
||
) |
This returns the generalized inverse of the given matrix.
This is the same as the inverse, except that if A is rank deficient we set the SVD value to 0.
This variation takes a boolean array indicating which parameters are to be held at zero.
A | Matrix to get inverse of |
Zero_unused | Boolean array indicating parameters that are unused and should be held to zero |
Rcond | We set singular values < Rcond * max singular value to 0. |
Definition at line 164 of file linear_algebra.cc.
blitz::Array< double, 2 > FullPhysics::multifit_covar | ( | const blitz::Array< double, 2 > & | J, |
double | Eps_rel = std::numeric_limits<double>::epsilon() |
||
) |
Definition at line 226 of file linear_algebra.cc.
bool FullPhysics::multifit_test_delta | ( | const blitz::Array< double, 1 > & | dx, |
const blitz::Array< double, 1 > & | x, | ||
double | Eps_abs = std::numeric_limits<double>::epsilon() , |
||
double | Eps_rel = std::numeric_limits<double>::epsilon() |
||
) |
Definition at line 239 of file linear_algebra.cc.
bool FullPhysics::multifit_test_gradient | ( | const blitz::Array< double, 1 > & | g, |
double | Eps_abs = std::numeric_limits<double>::epsilon() |
||
) |
Definition at line 234 of file linear_algebra.cc.
Array< double, 1 > FullPhysics::solve_least_squares | ( | const blitz::Array< double, 2 > & | A, |
const blitz::Array< double, 1 > & | B, | ||
double | Rcond = 1e-12 |
||
) |
This solves the least squares system A*x = b, returning the minimum norm solution.
It is perfectly ok for A to be rank deficient.
A | a M x N matrix, which may be rank deficient. |
B | the right hand side, a M vector |
Rcond | a singular value <= Rcond * max singular value is treated as 0. You can pass Rcond < 0 to use machine precision if desired. |
Definition at line 27 of file linear_algebra.cc.
Array< double, 1 > FullPhysics::solve_least_squares_qr | ( | const blitz::Array< double, 2 > & | A, |
const blitz::Array< double, 1 > & | B | ||
) |
This finds the least squares solution to the overdetermined system A x = b where the matrix A has more rows than columns.
The least squares solution minimizes the Euclidean norm of the residual, ||Ax - b||. The solution is determined using a QR decomposition of A.
A | a M x N matrix, which may be rank deficient. |
B | the right hand side, a M vector |
Definition at line 66 of file linear_algebra.cc.
|
inline |
Definition at line 11 of file linear_algebra.h.
void FullPhysics::svd | ( | const blitz::Array< double, 2 > & | A, |
blitz::Array< double, 1 > & | S, | ||
blitz::Array< double, 2 > & | U, | ||
blitz::Array< double, 2 > & | VT | ||
) |
Compute the SVD decomposition of a matrix A = U * SIGMA * V^T.
A | Matrix to get decomposition for. |
S | On exit, Singular values |
U | On exit, U matrix. Note that this is the "thin" version of U, this is M x N rather than M x M. The "full" version would just have extra columns of zeros. |
VT | On exit, V^T matrix |
Definition at line 103 of file linear_algebra.cc.
blitz::Array<T, D> FullPhysics::to_c_order | ( | const blitz::Array< T, D > & | In | ) |
Ensure that a given blitz::Array is contiguous, not reversed, and in C RowMajorArray format.
If it already is, then we just return the array. Otherwise, we create a copy of this in C order.
Definition at line 44 of file linear_algebra.h.
blitz::Array<T, D> FullPhysics::to_c_order_const | ( | const blitz::Array< T, D > & | In | ) |
Ensure that a given blitz::Array is contiguous, not reversed, and in C RowMajorArray format.
If it already is, then we just return the array. Otherwise, we create a copy of this in C order.
Unfortunately, this removes the "const" from the calling argument, so we have a version of this where you explicitly say this ok.
Definition at line 70 of file linear_algebra.h.
blitz::Array<T, D> FullPhysics::to_fortran | ( | const blitz::Array< T, D > & | In | ) |
Ensure that a given blitz::Array is contiguous, not reversed, and in fortran ColumnMajorArray format.
If it already is, then we just return the array. Otherwise, we create a copy of this suitable for passing to Fortran.
Definition at line 21 of file linear_algebra.h.
blitz::Array<T, D> FullPhysics::to_fortran_const | ( | const blitz::Array< T, D > & | In | ) |
Ensure that a given blitz::Array is contiguous, not reversed, and in fortran ColumnMajorArray format.
If it already is, then we just return the array. Otherwise, we create a copy of this suitable for passing to Fortran.
Unfortunately, this removes the "const" from the calling argument, so we have a version of this where you explicitly say this ok.
Definition at line 156 of file linear_algebra.h.