iterative-solver 0.0
gram_schmidt.h
1#ifndef LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_SUBSPACE_GRAM_SCHMIDT_H
2#define LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_SUBSPACE_GRAM_SCHMIDT_H
3#include <molpro/linalg/array/ArrayHandler.h>
4#include <molpro/linalg/itsolv/subspace/Matrix.h>
5#include <molpro/linalg/itsolv/wrap.h>
6
37template <typename T>
38std::vector<T> gram_schmidt(const Matrix<T>& s, Matrix<T>& l, double norm_thresh = 1.0e-14) {
39 assert(s.rows() == s.cols());
40 auto n = s.rows();
41 l.fill(0);
42 l.resize({n, n});
43 auto norm = std::vector<double>(n, 0);
44 auto w = std::vector<double>{};
45 for (size_t i = 0; i < n; ++i) {
46 w.assign(i, 0.);
47 for (size_t j = 0; j < i; ++j) {
48 for (size_t k = 0; k <= j; ++k) {
49 w[j] += s(i, k) * l(j, k);
50 }
51 }
52 for (size_t j = 0; j < i; ++j) {
53 if (norm[j] > norm_thresh) {
54 for (size_t k = 0; k <= j; ++k) {
55 l(i, k) -= w[j] / norm[j] * l(j, k);
56 }
57 }
58 }
59 l(i, i) = 1.;
60 for (size_t j = 0; j <= i; ++j) {
61 for (size_t k = 0; k < j; ++k) {
62 norm[i] += 2 * l(i, j) * l(i, k) * s(j, k);
63 }
64 norm[i] += l(i, j) * l(i, j) * s(j, j);
65 }
66 }
67 std::transform(begin(norm), end(norm), begin(norm), [](auto el) { return std::sqrt(std::abs(el)); });
68 return norm;
69}
70
71// FIXME This is implicitly constructed in Gram Schmidt and we can make it an opitonal return variable
97template <typename value_type, typename value_type_abs>
99 const Matrix<value_type>& lin_trans,
100 const std::vector<value_type_abs>& norm) {
101 const auto nrows = lin_trans.rows();
102 const auto ncols = lin_trans.cols();
103 auto t = Matrix<value_type>({nrows, ncols});
104 for (size_t i = 0; i < nrows; ++i) {
105 for (size_t j = 0; j < i; ++j) {
106 for (size_t k = 0; k <= j; ++k) {
107 t(i, j) -= lin_trans(j, k) * overlap(i, k) / std::pow(norm[j], 2);
108 }
109 }
110 }
111 return t;
112}
113
127template <class R, typename value_type_abs>
128auto modified_gram_schmidt(VecRef<R>& params, array::ArrayHandler<R, R>& handler, value_type_abs null_thresh) {
129 auto null_param_indices = std::vector<size_t>{};
130 const size_t n = params.size();
131 for (size_t i = 0; i < n; ++i) {
132 auto norm = handler.dot(params[i], params[i]);
133 norm = std::sqrt(std::abs(norm));
134 if (norm > null_thresh) {
135 handler.scal(1. / norm, params[i]);
136 for (size_t j = i + 1; j < n; ++j) {
137 auto ov = handler.dot(params[i], params[j]);
138 handler.axpy(-ov, params[i], params[j]);
139 }
140 } else {
141 null_param_indices.emplace_back(i);
142 }
143 }
144 return null_param_indices;
145}
146
147} // namespace molpro::linalg::itsolv::subspace::util
148
149#endif // LINEARALGEBRA_SRC_MOLPRO_LINALG_ITSOLV_SUBSPACE_GRAM_SCHMIDT_H
Enhances various operations between pairs of arrays and allows dynamic code injection with uniform in...
Definition: ArrayHandler.h:162
virtual value_type dot(const AL &x, const AR &y)=0
virtual void scal(value_type alpha, AL &x)=0
virtual void axpy(value_type alpha, const AR &x, AL &y)=0
Matrix container that allows simple data access, slicing, copying and resizing without loosing data.
Definition: Matrix.h:28
void fill(T value)
Sets all elements of matrix to value.
Definition: Matrix.h:89
void resize(const coord_type &dims)
Resize the matrix. The old data is preserved and any new rows/cols are zeroed.
Definition: Matrix.h:126
index_type rows() const
Definition: Matrix.h:165
index_type cols() const
Definition: Matrix.h:166
Definition: gram_schmidt.h:7
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.
Definition: util.h:49
std::vector< T > gram_schmidt(const Matrix< T > &s, Matrix< T > &l, double norm_thresh=1.0e-14)
Performs Gram-Schmidt orthogonalisation without normalisation.
Definition: gram_schmidt.h:38
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.
Definition: gram_schmidt.h:98
auto modified_gram_schmidt(VecRef< R > &params, array::ArrayHandler< R, R > &handler, value_type_abs null_thresh)
Apply modified Gram-Schmidt procedure to orthonormalise parameters.
Definition: gram_schmidt.h:128
std::vector< std::reference_wrapper< A > > VecRef
Definition: wrap.h:11