givensRotations.hh

Go to the documentation of this file.
1 
7 #ifndef givensRotations_hh
8 #define givensRotations_hh
9 
10 #include "sparseqr/sparseqr.hh"
11 #include "operator/matrix.hh"
12 
13 namespace sparseqr {
14 
15  // ******************************************************* GivensRotations **
16 
31  template<typename F>
32  class GivensRotations : public concepts::Operator<F> {
33  public:
34  template<class G>
36  const Qmatrix& q, bool transpose)
37  : concepts::Operator<F>(space.dim(), space.dim())
38  , q_(q), transpose_(transpose) {}
39  GivensRotations(uint dim, const Qmatrix& q, bool transpose)
40  : concepts::Operator<F>(dim, dim), q_(q), transpose_(transpose) {}
41  virtual void operator()(const concepts::Function<F>& fncY,
42  concepts::Function<F>& fncX);
46  void multiply(const concepts::Matrix<F>& fact,
47  concepts::Matrix<F>& dest) const;
52  concepts::Matrix<F>& dest) const;
53  protected:
54  virtual std::ostream& info(std::ostream& os) const;
55  private:
57  const Qmatrix& q_;
58  const bool transpose_;
59  };
60 
61  template<typename F>
63  (const concepts::Function<F>& fncY, concepts::Function<F>& fncX) {
64  conceptsAssert(fncY.dim() == this->dimY(), concepts::Assertion());
65  conceptsAssert(fncX.dim() == this->dimX(), concepts::Assertion());
66  fncX = fncY;
67  const J* j = q_.qend;
68  if (transpose_) j = q_.qbegin;
69  while (j != 0) {
70  conceptsAssert((uint)j->row1 < this->dimY(), concepts::Assertion());
71  conceptsAssert((uint)j->row2 < this->dimY(), concepts::Assertion());
72  F t = fncX(j->row1);
73  if (transpose_) {
74  fncX(j->row1) = j->c*fncX(j->row1) + j->s*fncX(j->row2);
75  fncX(j->row2) = j->c*fncX(j->row2) - j->s*t;
76  } else {
77  fncX(j->row1) = j->c*fncX(j->row1) - j->s*fncX(j->row2);
78  fncX(j->row2) = j->c*fncX(j->row2) + j->s*t;
79  }
80  if (transpose_) j = j->next;
81  else j = j->previous;
82  }
83  }
84 
85  /* The following routine does the same for a matrix as operator() above
86  does for a vector.
87 
88  It has to compute \c dest = \c fact times Given's rotations whereas
89  operator() computes \c fncX = Given's rotations times \c fncY.
90  \c dest^T = Given's^T times \c fact^T => each column of \c fact^T is
91  treated the same way as \c fncY. A column of \c fact^T is a row
92  of \c fact.
93 
94  The precondition of the routine is \c fact = \c dest, therefore, we
95  can work on \c dest only.
96  */
97  template<typename F>
99  (const concepts::Matrix<F>& fact, concepts::Matrix<F>& dest) const {
100  // check sizes of the matrices
101  conceptsAssert(fact.dimX() == dest.dimX(), concepts::Assertion());
102  conceptsAssert(fact.dimY() == dest.dimY(), concepts::Assertion());
103  conceptsAssert(fact.dimY() == this->dimX(), concepts::Assertion());
104  // loop over rows of dest
105  for (uint i = 0; i < dest.dimX(); ++i) {
106  const J* j = q_.qbegin;
107  if (transpose_) j = q_.qend;
108  while (j != 0) { // loop over Given's rotations
109  conceptsAssert((uint)j->row1 < this->dimX(), concepts::Assertion());
110  conceptsAssert((uint)j->row2 < this->dimX(), concepts::Assertion());
111  F t = dest(i, j->row1);
112  if (transpose_) {
113  dest(i, j->row1) = j->c*dest(i, j->row1) - j->s*dest(i, j->row2);
114  dest(i, j->row2) = j->c*dest(i, j->row2) + j->s*t;
115  } else {
116  dest(i, j->row1) = j->c*dest(i, j->row1) + j->s*dest(i, j->row2);
117  dest(i, j->row2) = j->c*dest(i, j->row2) - j->s*t;
118  }
119  if (transpose_) j = j->previous;
120  else j = j->next;
121  } // while j
122  } // for i
123  }
124 
125  template<typename F>
127  (const concepts::Matrix<F>& fact, concepts::Matrix<F>& dest) const {
128  // check sizes of the matrices
129  conceptsAssert(fact.dimX() == dest.dimX(), concepts::Assertion());
130  conceptsAssert(fact.dimY() == dest.dimY(), concepts::Assertion());
131  conceptsAssert(fact.dimX() == this->dimX(), concepts::Assertion());
132  // loop over columns of dest
133  for (uint i = 0; i < dest.dimY(); ++i) {
134  const J* j = q_.qend;
135  if (transpose_) j = q_.qbegin;
136  while (j != 0) {
137  conceptsAssert((uint)j->row1 < this->dimX(), concepts::Assertion());
138  conceptsAssert((uint)j->row2 < this->dimX(), concepts::Assertion());
139  F t = dest(j->row1, i);
140  if (transpose_) {
141  dest(j->row1, i) = j->c*dest(j->row1, i) + j->s*dest(j->row2, i);
142  dest(j->row2, i) = j->c*dest(j->row2, i) - j->s*t;
143  } else {
144  dest(j->row1, i) = j->c*dest(j->row1, i) - j->s*dest(j->row2, i);
145  dest(j->row2, i) = j->c*dest(j->row2, i) + j->s*t;
146  }
147  if (transpose_) j = j->next;
148  else j = j->previous;
149  }
150  }
151  }
152 
153  template<typename F>
154  std::ostream& GivensRotations<F>::info(std::ostream& os) const {
155  return os << concepts::typeOf(*this)<<"(" << this->dimX() << ", " << q_ << ')';
156  }
157 
158 } // namespace sparseqr
159 
160 #endif // givensRotations_hh
Operator(uint dimX, uint dimY)
Definition: compositions.hh:51
J * previous
Definition: sparseqr.hh:84
double c
Definition: sparseqr.hh:83
Abstract class for a function.
Definition: basis.hh:21
virtual std::ostream & info(std::ostream &os) const
Given's rotations.
Definition: sparseqr.hh:80
GivensRotations(const concepts::Space< G > &space, const Qmatrix &q, bool transpose)
#define conceptsAssert(cond, exc)
Assert that a certain condition is fulfilled.
Definition: exceptions.hh:394
double s
Definition: sparseqr.hh:83
virtual void operator()(const concepts::Function< F > &fncY, concepts::Function< F > &fncX)
SparseQR solver.
Definition: driver.hh:19
Exception class for assertions.
Definition: exceptions.hh:258
Abstract class for an operator.
Definition: compositions.hh:31
void multiply(const concepts::Matrix< F > &fact, concepts::Matrix< F > &dest) const
Computes fact times Givens and writes the result to dest.
Abstract class for an operator.
Definition: ARPACK.hh:16
Mapping< F, dim > & transpose(Mapping< F, dim > &m)
Definition: operations.hh:54
void multiplyFirst(const concepts::Matrix< F > &fact, concepts::Matrix< F > &dest) const
Computes Givens times fact and writes the result to dest.
std::string typeOf(const T &t)
Return the typeid name of a class object.
Definition: output.hh:43
const Qmatrix & q_
Given's rotations.
Q matrix of the QR factorization.
Definition: sparseqr.hh:96
GivensRotations(uint dim, const Qmatrix &q, bool transpose)
Basic namespace for Concepts-2.
Definition: pml_formula.h:16
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich