Functions | |
| template<typename Matrix > | |
| void | assert_left_cols_to_skip (Matrix &in_output, Eigen::Index left_cols_to_skip) |
| template<typename Matrix > | |
| void | GS_orthogonalisation (Matrix &in_output, Eigen::Index left_cols_to_skip=0) |
| template<typename Matrix > | |
| void | JensWehner_orthogonalisation (Matrix &in_output, Eigen::Index left_cols_to_skip=0) |
| template<typename Matrix > | |
| void | MGS_orthogonalisation (Matrix &in_output, Eigen::Index left_cols_to_skip=0) |
| template<typename Matrix > | |
| void | QR_orthogonalisation (Matrix &in_output) |
| template<typename Matrix > | |
| void | subspace_orthogonalisation (Matrix &in_output, Eigen::Index left_cols_to_skip) |
| template<typename Matrix > | |
| Eigen::Index | treat_first_col (Matrix &in_output, Eigen::Index left_cols_to_skip) |
| template<typename Matrix > | |
| void | twice_is_enough_orthogonalisation (Matrix &in_output, Eigen::Index left_cols_to_skip=0) |
| void Spectra::assert_left_cols_to_skip | ( | Matrix & | in_output, |
| Eigen::Index | left_cols_to_skip | ||
| ) |
Check if the number of columns to skip is larger than 0 but smaller than the total number of columns of the matrix
| in_output | Matrix to be orthogonalized |
| left_cols_to_skip | Number of left columns to be left untouched |
Definition at line 21 of file Orthogonalization.h.
| void Spectra::GS_orthogonalisation | ( | Matrix & | in_output, |
| Eigen::Index | left_cols_to_skip = 0 |
||
| ) |
Orthogonalize the in_output matrix using a Gram Schmidt process
| in_output | matrix to be orthogonalized |
| left_cols_to_skip | Number of left columns to be left untouched |
Definition at line 80 of file Orthogonalization.h.
| void Spectra::JensWehner_orthogonalisation | ( | Matrix & | in_output, |
| Eigen::Index | left_cols_to_skip = 0 |
||
| ) |
Orthogonalize the in_output matrix using a Jens process The subspace spanned by right columns are first orthogonalized agains the left columns, and then a QR decomposition is applied on the right columns to make them orthogonalized agains each other
| in_output | Matrix to be orthogonalized |
| left_cols_to_skip | Number of left columns to be left untouched |
Definition at line 119 of file Orthogonalization.h.
| void Spectra::MGS_orthogonalisation | ( | Matrix & | in_output, |
| Eigen::Index | left_cols_to_skip = 0 |
||
| ) |
Orthogonalize the in_output matrix using a modified Gram Schmidt process
| in_output | matrix to be orthogonalized |
| left_cols_to_skip | Number of left columns to be left untouched |
Definition at line 61 of file Orthogonalization.h.
| void Spectra::QR_orthogonalisation | ( | Matrix & | in_output | ) |
Orthogonalize the in_output matrix using a QR decomposition
| in_output | Matrix to be orthogonalized |
Definition at line 46 of file Orthogonalization.h.
| void Spectra::subspace_orthogonalisation | ( | Matrix & | in_output, |
| Eigen::Index | left_cols_to_skip | ||
| ) |
Orthogonalize the subspace spanned by right columns of in_output against the subspace spanned by left columns It assumes that the left columns are already orthogonal and normalized, and it does not orthogonalize the left columns against each other
| in_output | Matrix to be orthogonalized |
| left_cols_to_skip | Number of left columns to be left untouched |
Definition at line 99 of file Orthogonalization.h.
| Eigen::Index Spectra::treat_first_col | ( | Matrix & | in_output, |
| Eigen::Index | left_cols_to_skip | ||
| ) |
If the the number of columns to skip is null, normalize the first column and set left_cols_to_skip=1
| in_output | Matrix to be orthogonalized |
| left_cols_to_skip | Number of left columns to be left untouched |
Definition at line 33 of file Orthogonalization.h.
| void Spectra::twice_is_enough_orthogonalisation | ( | Matrix & | in_output, |
| Eigen::Index | left_cols_to_skip = 0 |
||
| ) |
Orthogonalize the in_output matrix using a twice-is-enough Jens process
| in_output | Matrix to be orthogonalized |
| left_cols_to_skip | Number of left columns to be left untouched |
Definition at line 133 of file Orthogonalization.h.