Class documentation of Concepts

Loading...
Searching...
No Matches
driver.hh
Go to the documentation of this file.
1
6#ifndef sparseqrDriver_hh
7#define sparseqrDriver_hh
8
9#include <memory>
10#include "basics/typedefs.hh"
12#include "operator/sparseMatrix.hh"
14#include "sparseqr/sparseqr.hh"
16
17using concepts::Real;
18
19namespace sparseqr {
20
21 // **************************************************************** Driver **
22
34 class Driver : concepts::Operator<Real> {
35 public:
42 : concepts::Operator<Real>(A.dimY(), A.dimX())
43 , rank_(0), A_(A), computed_(false)
44 , Prow_(nullptr), Pcol_(nullptr), Q_(nullptr), Qt_(nullptr), restr_(nullptr), ext_(nullptr)
45 , prow_(0), pcol_(0), qr_(nullptr) {}
46 virtual ~Driver();
60 inline int rank();
61
62 virtual void operator()(const concepts::Function<Real>& fncY,
64 void operator()(const concepts::Vector<Real>& fncY,
66 protected:
67 virtual std::ostream& info(std::ostream& os) const;
68 private:
69 void compute_();
70 int rank_;
72 bool computed_;
73 std::unique_ptr<concepts::Permutation<Real> > Prow_, Pcol_;
74 std::unique_ptr<sparseqr::GivensRotations<Real> > Q_, Qt_;
75 std::unique_ptr<concepts::TrivExtendRestrict<Real> > restr_, ext_;
76 concepts::Array<int> prow_, pcol_;
77 std::unique_ptr<QR> qr_;
79 };
80
82 if (!computed_) compute_();
83 return rank_;
84 }
85
86} // namespace sparseqr
87
88#endif // sparseqrDriver_hh
virtual const uint dimY() const
virtual const uint dimX() const
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
Driver(const concepts::SparseMatrix< Real > &A)
Definition driver.hh:41
sparseqr::GivensRotations< Real > * Qt()
Returns .
concepts::TrivExtendRestrict< Real > * restriction()
Returns the restriction from n to the last n-r degrees of freedom.
concepts::TrivExtendRestrict< Real > * extension()
Returns the extension from the last n-r to n degrees of freedom.
sparseqr::GivensRotations< Real > * Q()
Returns .
int rank()
Returns the rank of the matrix.
Definition driver.hh:81
concepts::Permutation< Real > * Pcol()
Returns .
concepts::Permutation< Real > * Prow()
Returns .
double Real
Definition typedefs.hh:39