19 #include <deal.II/base/exceptions.h> 20 #include <deal.II/base/std_cxx14/memory.h> 22 #include <deal.II/lac/lapack_full_matrix.h> 23 #include <deal.II/lac/lapack_templates.h> 24 #include <deal.II/lac/utilities.h> 26 #include <boost/signals2/signal.hpp> 28 DEAL_II_NAMESPACE_OPEN
41 template <
typename VectorType>
47 typedef typename VectorType::value_type
Number;
94 solve(Vector<Number> & x,
95 const Vector<Number> &y,
135 boost::signals2::connection
137 const std::function<
void(
const unsigned int i,
138 const unsigned int j,
139 const std::array<Number, 3> &csr)> &slot);
158 std::vector<std::unique_ptr<VectorType>>
columns;
174 boost::signals2::signal<void(
const unsigned int i,
175 const unsigned int j,
176 const std::array<Number, 3> &)>
237 template <
typename VectorType>
244 typedef typename VectorType::value_type
Number;
254 virtual ~QR() =
default;
346 template <
typename VectorType>
353 typedef typename VectorType::value_type
Number;
404 boost::signals2::connection
406 const std::function<
bool(
const Vector<Number> &u,
408 const Number & col_norm_sqr)> &slot);
425 boost::signals2::signal<bool(const Vector<Number> &u,
427 const Number & col_norm_sqr)>
436 template <
typename VectorType>
445 template <
typename VectorType>
454 template <
typename VectorType>
463 template <
typename VectorType>
466 const Vector<Number> &y,
469 Assert(x.size() == this->current_size,
471 Assert(y.size() == this->current_size,
478 const int lda = this->current_size;
479 const int ldb = this->current_size;
480 const int N = this->current_size;
497 template <
typename VectorType>
500 const Vector<Number> &x)
const 502 Assert(x.size() == this->current_size,
506 for (
unsigned int j = 0; j < this->current_size; ++j)
507 y.add(x[j], *this->columns[j]);
512 template <
typename VectorType>
515 const VectorType &x)
const 517 Assert(y.size() == this->current_size,
520 for (
unsigned int j = 0; j < this->current_size; ++j)
521 y[j] = (*this->columns[j]) * x;
526 template <
class VectorType>
527 boost::signals2::connection
529 const std::function<
void(
const unsigned int i,
530 const unsigned int j,
531 const std::array<Number, 3> &)> &slot)
533 return givens_signal.connect(slot);
538 template <
class VectorType>
539 boost::signals2::connection
541 const std::function<
bool(
const Vector<Number> &u,
543 const Number & col_norm_sqr)> &slot)
545 return column_signal.connect(slot);
550 template <
typename VectorType>
557 template <
typename VectorType>
561 if (this->current_size == 0)
563 this->R.grow_or_shrink(this->current_size + 1);
564 this->columns.push_back(std_cxx14::make_unique<VectorType>(column));
565 this->R(0, 0) = column.l2_norm();
566 ++this->current_size;
571 Vector<Number> u(this->current_size);
572 this->multiply_with_AT(u, column);
575 const int lda = this->current_size;
576 const int ldb = this->current_size;
577 const int N = this->current_size;
581 "U",
"T",
"N", &N, &n_rhs, &this->R(0, 0), &lda, &u(0), &ldb, &info);
585 const Number column_norm_sqr = column.norm_sqr();
586 const Number rho2 = column_norm_sqr - u.norm_sqr();
587 const bool linearly_independent =
588 column_signal.empty() ? rho2 > 0 :
589 column_signal(u, rho2, column_norm_sqr).get();
592 if (!linearly_independent)
597 this->columns.push_back(std_cxx14::make_unique<VectorType>(column));
598 this->R.grow_or_shrink(this->current_size + 1);
599 this->R(this->current_size, this->current_size) = std::sqrt(rho2);
600 for (
unsigned int i = 0; i < this->current_size; ++i)
601 this->R(i, this->current_size) = u(i);
603 this->current_size++;
611 template <
typename VectorType>
614 const unsigned int k)
618 const std::array<Number, 3> csr =
619 ::Utilities::LinearAlgebra::givens_rotation<Number>(this->R(i, k),
623 this->R(i, k) = csr[2];
626 for (
unsigned int j = 0; j < this->R.n(); ++j)
629 const Number t = this->R(i, j);
630 this->R(i, j) = csr[0] * this->R(i, j) + csr[1] * this->R(k, j);
631 this->R(k, j) = -csr[1] * t + csr[0] * this->R(k, j);
634 if (!this->givens_signal.empty())
635 this->givens_signal(i, k, csr);
640 template <
typename VectorType>
646 for (
unsigned int j = k + 1; j < this->R.n(); ++j)
648 const unsigned int i = j - 1;
649 apply_givens_rotation(i, j);
653 --this->current_size;
654 this->R.remove_row_and_column(this->current_size, k);
657 this->columns.erase(this->columns.begin() + k);
662 template <
typename VectorType>
665 const Vector<Number> &x)
const 669 Vector<Number> x1 = x;
671 multiply_with_A(y, x1);
676 template <
typename VectorType>
679 const VectorType &x)
const 684 multiply_with_AT(y, x);
690 template <
typename VectorType>
693 const Vector<Number> &x)
const 700 template <
typename VectorType>
703 const VectorType &x)
const 710 template <
typename VectorType>
717 template <
typename VectorType>
722 this->R.grow_or_shrink(this->current_size + 1);
723 this->columns.push_back(std_cxx14::make_unique<VectorType>(column));
727 auto &last_col = *this->columns.back();
728 for (
unsigned int i = 0; i < this->current_size; ++i)
730 const auto &i_col = *this->columns[i];
731 this->R(i, this->current_size) = i_col * last_col;
732 last_col.add(-this->R(i, this->current_size), i_col);
735 this->R(this->current_size, this->current_size) = last_col.l2_norm();
737 Assert(this->R(this->current_size, this->current_size) > 0.,
739 last_col *= 1. / this->R(this->current_size, this->current_size);
741 ++this->current_size;
747 template <
typename VectorType>
750 const unsigned int k)
754 const std::array<Number, 3> csr =
755 ::Utilities::LinearAlgebra::givens_rotation<Number>(this->R(i, k),
759 this->R(i, k) = csr[2];
762 for (
unsigned int j = 0; j < this->R.n(); ++j)
765 const Number t = this->R(i, j);
766 this->R(i, j) = csr[0] * this->R(i, j) + csr[1] * this->R(k, j);
767 this->R(k, j) = -csr[1] * t + csr[0] * this->R(k, j);
772 auto &col_i = *this->columns[i];
773 auto &col_k = *this->columns[k];
776 col_i.sadd(csr[0], csr[1], col_k);
777 col_k.sadd(csr[0], -csr[1], tmp);
779 if (!this->givens_signal.empty())
780 this->givens_signal(i, k, csr);
785 template <
typename VectorType>
790 Assert(this->current_size > 0,
791 ExcMessage(
"Can not remove a column if QR is empty"));
808 for (
unsigned int j = k + 1; j < this->R.n(); ++j)
810 const unsigned int i = j - 1;
811 apply_givens_rotation(i, j);
817 const unsigned int size_minus_1 = this->columns.size() - 1;
818 this->columns.erase(this->columns.begin() + size_minus_1);
821 --this->current_size;
822 this->R.remove_row_and_column(this->current_size, k);
827 template <
typename VectorType>
836 template <
typename VectorType>
845 template <
typename VectorType>
849 Vector<Number> x1 = x;
850 const int N = this->current_size;
853 trmv(
"U",
"N",
"N", &N, &this->R(0, 0), &lda, &x1[0], &incx);
855 multiply_with_Q(y, x1);
860 template <
typename VectorType>
864 multiply_with_QT(y, x);
866 const int N = this->current_size;
869 trmv(
"U",
"T",
"N", &N, &this->R(0, 0), &lda, &y[0], &incx);
874 DEAL_II_NAMESPACE_CLOSE
unsigned int size() const
VectorType::value_type Number
std::vector< std::unique_ptr< VectorType > > columns
Matrix is upper triangular.
void multiply_with_colsT(Vector< Number > &y, const VectorType &x) const
virtual bool append_column(const VectorType &column)
virtual ~BaseQR()=default
virtual void multiply_with_A(VectorType &y, const Vector< Number > &x) const =0
unsigned int current_size
virtual void multiply_with_AT(Vector< Number > &y, const VectorType &x) const
virtual void multiply_with_Q(VectorType &y, const Vector< Number > &x) const =0
virtual bool append_column(const VectorType &column)=0
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual void multiply_with_Q(VectorType &y, const Vector< Number > &x) const
static ::ExceptionBase & ExcDivideByZero()
LAPACKFullMatrix< Number > R
virtual void multiply_with_A(VectorType &y, const Vector< Number > &x) const
boost::signals2::signal< bool(const Vector< Number > &u, const Number &rho, const Number &col_norm_sqr)> column_signal
boost::signals2::connection connect_append_column_slot(const std::function< bool(const Vector< Number > &u, const Number &rho2, const Number &col_norm_sqr)> &slot)
static ::ExceptionBase & ExcMessage(std::string arg1)
void solve(Vector< Number > &x, const Vector< Number > &y, const bool transpose=false) const
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)
virtual ~ImplicitQR()=default
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
VectorType::value_type Number
virtual void multiply_with_AT(Vector< Number > &y, const VectorType &x) const =0
virtual void remove_column(const unsigned int k=0)
boost::signals2::signal< void(const unsigned int i, const unsigned int j, const std::array< Number, 3 > &)> givens_signal
virtual void multiply_with_A(VectorType &y, const Vector< Number > &x) const
virtual void multiply_with_AT(Vector< Number > &y, const VectorType &x) const
const LAPACKFullMatrix< Number > & get_R() const
virtual void multiply_with_Q(VectorType &y, const Vector< Number > &x) const
SymmetricTensor< rank_, dim, Number > transpose(const SymmetricTensor< rank_, dim, Number > &t)
void apply_givens_rotation(const unsigned int i, const unsigned int k)
virtual void multiply_with_QT(Vector< Number > &y, const VectorType &x) const
virtual bool append_column(const VectorType &column)
virtual void remove_column(const unsigned int k=0)
virtual void multiply_with_QT(Vector< Number > &y, const VectorType &x) const
virtual void remove_column(const unsigned int k=0)=0
VectorType::value_type Number
void apply_givens_rotation(const unsigned int i, const unsigned int k)
void set_property(const LAPACKSupport::Property property)
virtual void multiply_with_QT(Vector< Number > &y, const VectorType &x) const =0
void multiply_with_cols(VectorType &y, const Vector< Number > &x) const