26#include <boost/signals2/signal.hpp>
43template <
typename VectorType>
49 using Number =
typename VectorType::value_type;
137 boost::signals2::connection
139 const std::function<
void(
const unsigned int i,
140 const unsigned int j,
141 const std::array<Number, 3> &csr)> &slot);
160 std::vector<std::unique_ptr<VectorType>>
columns;
176 boost::signals2::signal<void(
const unsigned int i,
177 const unsigned int j,
178 const std::array<Number, 3> &)>
239template <
typename VectorType>
246 using Number =
typename VectorType::value_type;
348template <
typename VectorType>
355 using Number =
typename VectorType::value_type;
406 boost::signals2::connection
410 const Number & col_norm_sqr)> &slot);
429 const Number & col_norm_sqr)>
438 namespace QRImplementation
445 template <
typename Number>
456 template <
typename Number>
473template <
typename VectorType>
482template <
typename VectorType>
491template <
typename VectorType>
500template <
typename VectorType>
515 const int lda = this->current_size;
516 const int ldb = this->current_size;
517 const int N = this->current_size;
534template <
typename VectorType>
543 for (
unsigned int j = 0; j < this->current_size; ++j)
544 y.add(x[j], *this->columns[j]);
549template <
typename VectorType>
552 const VectorType &x)
const
557 for (
unsigned int j = 0; j < this->current_size; ++j)
558 y[j] = (*this->columns[j]) * x;
563template <
class VectorType>
564boost::signals2::connection
566 const std::function<
void(
const unsigned int i,
567 const unsigned int j,
568 const std::array<Number, 3> &)> &slot)
570 return givens_signal.connect(slot);
575template <
class VectorType>
576boost::signals2::connection
580 const Number & col_norm_sqr)> &slot)
582 return column_signal.connect(slot);
587template <
typename VectorType>
594template <
typename VectorType>
598 if (this->current_size == 0)
600 this->R.grow_or_shrink(this->current_size + 1);
601 this->columns.push_back(std::make_unique<VectorType>(column));
602 this->R(0, 0) = column.l2_norm();
603 ++this->current_size;
609 this->multiply_with_AT(u, column);
612 const int lda = this->current_size;
613 const int ldb = this->current_size;
614 const int N = this->current_size;
618 'U',
'T',
'N', N, n_rhs, &this->R(0, 0), lda, &u(0), ldb, &info);
622 const Number column_norm_sqr = column.norm_sqr();
623 const Number rho2 = column_norm_sqr - u.
norm_sqr();
624 const bool linearly_independent =
625 column_signal.empty() ? rho2 > 0 :
626 column_signal(u, rho2, column_norm_sqr).get();
629 if (!linearly_independent)
634 this->columns.push_back(std::make_unique<VectorType>(column));
635 this->R.grow_or_shrink(this->current_size + 1);
636 this->R(this->current_size, this->current_size) =
std::sqrt(rho2);
637 for (
unsigned int i = 0; i < this->current_size; ++i)
638 this->R(i, this->current_size) = u(i);
640 this->current_size++;
648template <
typename VectorType>
651 const unsigned int k)
655 const std::array<Number, 3> csr =
656 ::Utilities::LinearAlgebra::givens_rotation<Number>(this->R(i, k),
660 this->R(i, k) = csr[2];
663 for (
unsigned int j = 0; j < this->R.n(); ++j)
666 const Number t = this->R(i, j);
667 this->R(i, j) = csr[0] * this->R(i, j) + csr[1] * this->R(k, j);
668 this->R(k, j) = -csr[1] * t + csr[0] * this->R(k, j);
671 if (!this->givens_signal.empty())
672 this->givens_signal(i, k, csr);
677template <
typename VectorType>
683 for (
unsigned int j = k + 1; j < this->R.n(); ++j)
685 const unsigned int i = j - 1;
686 apply_givens_rotation(i, j);
690 --this->current_size;
691 this->R.remove_row_and_column(this->current_size, k);
694 this->columns.erase(this->columns.begin() + k);
699template <
typename VectorType>
708 multiply_with_A(y, x1);
713template <
typename VectorType>
716 const VectorType &x)
const
721 multiply_with_AT(y, x);
727template <
typename VectorType>
737template <
typename VectorType>
740 const VectorType &x)
const
747template <
typename VectorType>
754template <
typename VectorType>
759 this->R.grow_or_shrink(this->current_size + 1);
760 this->columns.push_back(std::make_unique<VectorType>(column));
764 auto &last_col = *this->columns.back();
765 for (
unsigned int i = 0; i < this->current_size; ++i)
767 const auto &i_col = *this->columns[i];
768 this->R(i, this->current_size) = i_col * last_col;
769 last_col.add(-this->R(i, this->current_size), i_col);
772 this->R(this->current_size, this->current_size) = last_col.l2_norm();
774 Assert(this->R(this->current_size, this->current_size) > 0.,
776 last_col *= 1. / this->R(this->current_size, this->current_size);
778 ++this->current_size;
784template <
typename VectorType>
787 const unsigned int k)
791 const std::array<Number, 3> csr =
792 ::Utilities::LinearAlgebra::givens_rotation<Number>(this->R(i, k),
796 this->R(i, k) = csr[2];
799 for (
unsigned int j = 0; j < this->R.n(); ++j)
802 const Number t = this->R(i, j);
803 this->R(i, j) = csr[0] * this->R(i, j) + csr[1] * this->R(k, j);
804 this->R(k, j) = -csr[1] * t + csr[0] * this->R(k, j);
809 auto &col_i = *this->columns[i];
810 auto &col_k = *this->columns[k];
813 col_i.sadd(csr[0], csr[1], col_k);
814 col_k.sadd(csr[0], -csr[1], tmp);
816 if (!this->givens_signal.empty())
817 this->givens_signal(i, k, csr);
822template <
typename VectorType>
827 Assert(this->current_size > 0,
828 ExcMessage(
"Can not remove a column if QR is empty"));
845 for (
unsigned int j = k + 1; j < this->R.n(); ++j)
847 const unsigned int i = j - 1;
848 apply_givens_rotation(i, j);
854 const unsigned int size_minus_1 = this->columns.size() - 1;
855 this->columns.erase(this->columns.begin() + size_minus_1);
858 --this->current_size;
859 this->R.remove_row_and_column(this->current_size, k);
864template <
typename VectorType>
873template <
typename VectorType>
882template <
typename VectorType>
887 const int N = this->current_size;
891 'U',
'N',
'N', N, &this->R(0, 0), lda, &x1[0], incx);
893 multiply_with_Q(y, x1);
898template <
typename VectorType>
902 multiply_with_QT(y, x);
904 const int N = this->current_size;
908 'U',
'T',
'N', N, &this->R(0, 0), lda, &y[0], incx);
boost::signals2::connection connect_givens_slot(const std::function< void(const unsigned int i, const unsigned int j, const std::array< Number, 3 > &csr)> &slot)
void multiply_with_colsT(Vector< Number > &y, const VectorType &x) const
unsigned int size() const
virtual void remove_column(const unsigned int k=0)=0
virtual ~BaseQR()=default
virtual void multiply_with_Q(VectorType &y, const Vector< Number > &x) const =0
virtual void multiply_with_A(VectorType &y, const Vector< Number > &x) const =0
std::vector< std::unique_ptr< VectorType > > columns
virtual void multiply_with_AT(Vector< Number > &y, const VectorType &x) const =0
void solve(Vector< Number > &x, const Vector< Number > &y, const bool transpose=false) const
LAPACKFullMatrix< Number > R
typename VectorType::value_type Number
virtual bool append_column(const VectorType &column)=0
boost::signals2::signal< void(const unsigned int i, const unsigned int j, const std::array< Number, 3 > &)> givens_signal
void multiply_with_cols(VectorType &y, const Vector< Number > &x) const
unsigned int current_size
const LAPACKFullMatrix< Number > & get_R() const
virtual void multiply_with_QT(Vector< Number > &y, const VectorType &x) const =0
virtual ~ImplicitQR()=default
virtual void multiply_with_QT(Vector< Number > &y, const VectorType &x) const
virtual void multiply_with_A(VectorType &y, const Vector< Number > &x) const
boost::signals2::connection connect_append_column_slot(const std::function< bool(const Vector< Number > &u, const Number &rho2, const Number &col_norm_sqr)> &slot)
typename VectorType::value_type Number
virtual bool append_column(const VectorType &column)
virtual void multiply_with_Q(VectorType &y, const Vector< Number > &x) const
virtual void remove_column(const unsigned int k=0)
boost::signals2::signal< bool(const Vector< Number > &u, const Number &rho, const Number &col_norm_sqr)> column_signal
void apply_givens_rotation(const unsigned int i, const unsigned int k)
virtual void multiply_with_AT(Vector< Number > &y, const VectorType &x) const
void set_property(const LAPACKSupport::Property property)
virtual void multiply_with_QT(Vector< Number > &y, const VectorType &x) const
typename VectorType::value_type Number
virtual bool append_column(const VectorType &column)
void apply_givens_rotation(const unsigned int i, const unsigned int k)
virtual void multiply_with_A(VectorType &y, const Vector< Number > &x) const
virtual void multiply_with_AT(Vector< Number > &y, const VectorType &x) const
virtual void remove_column(const unsigned int k=0)
virtual void multiply_with_Q(VectorType &y, const Vector< Number > &x) const
real_type norm_sqr() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcDivideByZero()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
@ upper_triangular
Matrix is upper triangular.
void call_trmv(const char uplo, const char trans, const char diag, const types::blas_int n, const Number *a, const types::blas_int lda, Number *x, const types::blas_int incx)
void call_trtrs(const char uplo, const char trans, const char diag, const types::blas_int n, const types::blas_int nrhs, const Number *a, const types::blas_int lda, Number *b, const types::blas_int ldb, types::blas_int *info)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)