Namespaces | |
namespace | detail |
Functions | |
template<typename T > | |
std::vector< T > | gram_schmidt (const Matrix< T > &s, Matrix< T > &l, double norm_thresh=1.0e-14) |
Performs Gram-Schmidt orthogonalisation without normalisation. More... | |
template<typename value_type , typename value_type_abs > | |
Matrix< value_type > | construct_lin_trans_in_orthogonal_set (const Matrix< value_type > &overlap, const Matrix< value_type > &lin_trans, const std::vector< value_type_abs > &norm) |
Construct Gram-Schmidt linear transformation in orthogonal vectors. More... | |
template<class R , typename value_type_abs > | |
auto | modified_gram_schmidt (VecRef< R > ¶ms, array::ArrayHandler< R, R > &handler, value_type_abs null_thresh) |
Apply modified Gram-Schmidt procedure to orthonormalise parameters. More... | |
template<class R , class Q , class Z , class W > | |
auto | overlap (const CVecRef< R > &left, const CVecRef< Q > &right, array::ArrayHandler< Z, W > &handler) -> std::enable_if_t< detail::Z_and_W_are_one_of_R_and_Q< R, Q, Z, W >, Matrix< double > > |
Calculates overlap matrix between left and right vectors. More... | |
template<class R > | |
Matrix< double > | overlap (const CVecRef< R > ¶ms, array::ArrayHandler< R, R > &handler) |
Calculates overlap matrix for a parameter set. Matrix is symmetric by construction. More... | |
template<typename T > | |
void | matrix_symmetrize (Matrix< T > &mat) |
template<typename T > | |
Matrix< T >::coord_type | max_element_index (const std::list< size_t > &rows, const std::list< size_t > &cols, const Matrix< T > &mat) |
Return maximum element in a matrix along specified rows and columns. More... | |
template<typename Slice > | |
std::vector< size_t > | eye_order (const Slice &mat) |
Returns order of rows in a matrix slice that brings it closest to identity. More... | |
Matrix< value_type > molpro::linalg::itsolv::subspace::util::construct_lin_trans_in_orthogonal_set | ( | const Matrix< value_type > & | overlap, |
const Matrix< value_type > & | lin_trans, | ||
const std::vector< value_type_abs > & | norm | ||
) |
Construct Gram-Schmidt linear transformation in orthogonal vectors.
Let {v} be the original vectors, and {u} their orthogonal set.
We can construct {u} from {v} using lower triangular linear transformation matrix L u_0 = v_0 u_1 = v_1 + L_{1,0} v_0 u_i = \sum_{j=0}^i L_{i,j} v_j
However, we can also construct u_i using v_i and previous {u}, using a lower triangular transformation matrix T, u_i = v_i - \sum_{j=0}^{i-1} <v_i,u_j> / <u_j, u_j> u_j = v_i + \sum_j=0^{i-1} T_{i,j} u_j
T_{i,j} = - <v_i,u_j> / <u_j, u_j> = - \sum_{k=0}^{j} L_{j,k} <v_i,v_k> / <u_j, u_j>
T_{i,j} = 0 if i >= j
overlap | Overlap matrix of original non-orthogonal vectors |
lin_trans | Gram-Schmidt linear transformation to an orthogonal set in terms of original vectors |
norm | norms of orthogonal vectors constructed from lin_trans @Returns linear transformation set |
std::vector< size_t > molpro::linalg::itsolv::subspace::util::eye_order | ( | const Slice & | mat | ) |
Returns order of rows in a matrix slice that brings it closest to identity.
Slice | Slice type |
mat | square matrix |
std::vector< T > molpro::linalg::itsolv::subspace::util::gram_schmidt | ( | const Matrix< T > & | s, |
Matrix< T > & | l, | ||
double | norm_thresh = 1.0e-14 |
||
) |
Performs Gram-Schmidt orthogonalisation without normalisation.
Any orthogonal vector with a norm less than threshold is not orthogonalised against.
For a vector set {v} with overlap matrix , generate a linear transformation L to a new orthogonal vector set
.
u_0 = v_0 u_1 = v_1 - <v_1, u_0> / <u_0, u_0> u_i = v_i -\sum_j=0^{i-1} <v_i, u_j> / <u_j, u_j> u_j = v_i -\sum_j=0^{i-1} <v_i, u_j> / <u_j, u_j> \sum_{k=0}^j L_{jk} v_k = v_i -\sum_j=0^{i-1} (\sum_{k=0}^{j} W_{ij} / N_{jj} L_{jk}) v_k
Let the overlap between u and v be,
W_{ij} = <v_i, u_j> = \sum_{k=0}^{j} S_{ik} L_{jk} N_{ii} = <u_i, u_i> = \sum_{k} L_{ik} <v_k, u_i> = \sum_{km} L_{ik} L_{im} S_{km}
L_{0,j} = \delta_{i,j} L_{i,k} = \delta_{i,k} - (1 + \delta_{i,k})(\sum_{j=0}^{k} W_{ij} / N_{jj} L_{jk})
s | overlap matrix |
l | row matrix with linear transformation into orthonormal basis. Upper triangular component is zero. |
norm_thresh | generated vector with norm less than threshold is not used for subsequent orthogonalisation |
void molpro::linalg::itsolv::subspace::util::matrix_symmetrize | ( | Matrix< T > & | mat | ) |
Matrix< T >::coord_type molpro::linalg::itsolv::subspace::util::max_element_index | ( | const std::list< size_t > & | rows, |
const std::list< size_t > & | cols, | ||
const Matrix< T > & | mat | ||
) |
Return maximum element in a matrix along specified rows and columns.
auto molpro::linalg::itsolv::subspace::util::modified_gram_schmidt | ( | VecRef< R > & | params, |
array::ArrayHandler< R, R > & | handler, | ||
value_type_abs | null_thresh | ||
) |
Apply modified Gram-Schmidt procedure to orthonormalise parameters.
New vectors with norm less than threshold are considered null and are not normalised. Their indices in params are returned.
R | |
value_type_abs |
params | |
handler | |
null_thresh |
auto molpro::linalg::itsolv::subspace::util::overlap | ( | const CVecRef< R > & | left, |
const CVecRef< Q > & | right, | ||
array::ArrayHandler< Z, W > & | handler | ||
) | -> std::enable_if_t<detail::Z_and_W_are_one_of_R_and_Q<R, Q, Z, W>, Matrix<double>> |
Calculates overlap matrix between left and right vectors.
Matrix< double > molpro::linalg::itsolv::subspace::util::overlap | ( | const CVecRef< R > & | params, |
array::ArrayHandler< R, R > & | handler | ||
) |
Calculates overlap matrix for a parameter set. Matrix is symmetric by construction.