![]() |
deal.II version 9.7.1
|
#include <deal.II/lac/qr.h>
A base class for thin QR implementations.
This class and classes derived from it are meant to build 





As a consequence, matrices which have the same number of rows as each vector (i.e. 

Public Member Functions | |
| virtual | ~BaseQR ()=default |
| virtual bool | append_column (const VectorType &column)=0 |
| virtual void | remove_column (const unsigned int k=0)=0 |
| unsigned int | size () const |
| const LAPACKFullMatrix< Number > & | get_R () const |
| void | solve (Vector< Number > &x, const Vector< Number > &y, const bool transpose=false) const |
| virtual void | multiply_with_Q (VectorType &y, const Vector< Number > &x) const =0 |
| virtual void | multiply_with_QT (Vector< Number > &y, const VectorType &x) const =0 |
| virtual void | multiply_with_A (VectorType &y, const Vector< Number > &x) const =0 |
| virtual void | multiply_with_AT (Vector< Number > &y, const VectorType &x) const =0 |
| 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) |
Protected Member Functions | |
| BaseQR () | |
| void | multiply_with_cols (VectorType &y, const Vector< Number > &x) const |
| void | multiply_with_colsT (Vector< Number > &y, const VectorType &x) const |
Protected Attributes | |
| std::vector< std::unique_ptr< VectorType > > | columns |
| LAPACKFullMatrix< Number > | R |
| unsigned int | current_size |
| boost::signals2::signal< void(const unsigned int i, const unsigned int j, const std::array< Number, 3 > &)> | givens_signal |
Private Types | |
| using | Number = typename VectorType::value_type |
|
private |
|
protected |
Default private constructor.
|
pure virtual |
Append column to the QR factorization. Returns true if the result is successful, i.e. the columns are linearly independent. Otherwise the column is rejected and the return value is false.
Implemented in ImplicitQR< VectorType >, and QR< VectorType >.
|
pure virtual |
Remove a column k and update QR factorization.
Implemented in ImplicitQR< VectorType >, and QR< VectorType >.
Return size of the subspace.
| const LAPACKFullMatrix< Number > & BaseQR< VectorType >::get_R | ( | ) | const |
Return the current upper triangular matrix R.
| void BaseQR< VectorType >::solve | ( | Vector< Number > & | x, |
| const Vector< Number > & | y, | ||
| const bool | transpose = false ) const |
Solve 
x and y should be consistent with the current size of the subspace. If transpose is true, 
|
pure virtual |
Set 

Implemented in ImplicitQR< VectorType >, and QR< VectorType >.
|
pure virtual |
Set 

Implemented in ImplicitQR< VectorType >, and QR< VectorType >.
|
pure virtual |
Set 

Implemented in ImplicitQR< VectorType >, and QR< VectorType >.
|
pure virtual |
Set 

Implemented in ImplicitQR< VectorType >, and QR< VectorType >.
| boost::signals2::connection BaseQR< VectorType >::connect_givens_slot | ( | const std::function< void(const unsigned int i, const unsigned int j, const std::array< Number, 3 > &csr)> & | slot | ) |
Connect a slot to retrieve a notification when the Givens rotations are performed.
The function takes two indices, i and j, describing the plane of rotation, and a triplet of numbers csr (cosine, sine and radius, see Utilities::LinearAlgebra::givens_rotation()) which represents the rotation matrix.
|
protected |
Compute 

|
protected |
Multiply with transpose columns stored in the object.
|
protected |
|
protected |