OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
QR.hpp
1#ifndef JAMA_QR_H
2#define JAMA_QR_H
3
4#include <cmath>
5
7
8using namespace OpenMD;
9using namespace std;
10
11namespace JAMA {
12
13 /**
14 <p>
15 Classical QR Decompisition:
16 for an m-by-n matrix A with m >= n, the QR decomposition is an m-by-n
17 orthogonal matrix Q and an n-by-n upper triangular matrix R so that
18 A = Q*R.
19 <P>
20 The QR decompostion always exists, even if the matrix does not have
21 full rank, so the constructor will never fail. The primary use of the
22 QR decomposition is in the least squares solution of nonsquare systems
23 of simultaneous linear equations. This will fail if isFullRank()
24 returns 0 (false).
25
26 <p>
27 The Q and R factors can be retrived via the getQ() and getR()
28 methods. Furthermore, a solve() method is provided to find the
29 least squares solution of Ax=b using the QR factors.
30
31 <p>
32 (Adapted from JAMA, a Java Matrix Library, developed by jointly
33 by the Mathworks and NIST; see http://math.nist.gov/javanumerics/jama).
34 */
35
36 template<class Real>
37 class QR {
38 /** Array for internal storage of decomposition.
39 @serial internal array storage.
40 */
41
43
44 /** Row and column dimensions.
45 @serial column dimension.
46 @serial row dimension.
47 */
48 int m, n;
49
50 /** Array for internal storage of diagonal of R.
51 @serial diagonal of R.
52 */
54
55 public:
56 /**
57 Create a QR factorization object for A.
58
59 @param A rectangular (m>=n) matrix.
60 */
61 QR(const DynamicRectMatrix<Real>& A) /* constructor */
62 {
64 m = A.getNRow();
65 n = A.getNCol();
66 Rdiag = DynamicVector<Real>(n);
67 int i = 0, j = 0, k = 0;
68
69 // Main loop.
70 for (k = 0; k < n; k++) {
71 // Compute 2-norm of k-th column without under/overflow.
72 Real nrm = 0;
73 for (i = k; i < m; i++) {
74 nrm = hypot(nrm, QR_(i, k));
75 }
76
77 if (nrm != 0.0) {
78 // Form k-th Householder vector.
79 if (QR_(k, k) < 0) { nrm = -nrm; }
80 for (i = k; i < m; i++) {
81 QR_(i, k) /= nrm;
82 }
83 QR_(k, k) += 1.0;
84
85 // Apply transformation to remaining columns.
86 for (j = k + 1; j < n; j++) {
87 Real s = 0.0;
88 for (i = k; i < m; i++) {
89 s += QR_(i, k) * QR_(i, j);
90 }
91 s = -s / QR_(k, k);
92 for (i = k; i < m; i++) {
93 QR_(i, j) += s * QR_(i, k);
94 }
95 }
96 }
97 Rdiag(k) = -nrm;
98 }
99 }
100
101 /**
102 Flag to denote the matrix is of full rank.
103
104 @return 1 if matrix is full rank, 0 otherwise.
105 */
106 int isFullRank() const {
107 for (int j = 0; j < n; j++) {
108 if (Rdiag(j) == 0) return 0;
109 }
110 return 1;
111 }
112
113 /**
114
115 Retreive the Householder vectors from QR factorization
116 @returns lower trapezoidal matrix whose columns define the reflections
117 */
120
121 /* note: H is completely filled in by algorithm, so
122 initializaiton of H is not necessary.
123 */
124 for (int i = 0; i < m; i++) {
125 for (int j = 0; j < n; j++) {
126 if (i >= j) {
127 H(i, j) = QR_(i, j);
128 } else {
129 H(i, j) = 0.0;
130 }
131 }
132 }
133 return H;
134 }
135
136 /** Return the upper triangular factor, R, of the QR factorization
137 @return R
138 */
141 for (int i = 0; i < n; i++) {
142 for (int j = 0; j < n; j++) {
143 if (i < j) {
144 R(i, j) = QR_(i, j);
145 } else if (i == j) {
146 R(i, j) = Rdiag(i);
147 } else {
148 R(i, j) = 0.0;
149 }
150 }
151 }
152 return R;
153 }
154
155 /**
156 Generate and return the (economy-sized) orthogonal factor
157 @return Q the (economy-sized) orthogonal factor (Q*R=A).
158 */
160 int i = 0, j = 0, k = 0;
161
163 for (k = n - 1; k >= 0; k--) {
164 for (i = 0; i < m; i++) {
165 Q(i, k) = 0.0;
166 }
167 Q(k, k) = 1.0;
168 for (j = k; j < n; j++) {
169 if (QR_(k, k) != 0) {
170 Real s = 0.0;
171 for (i = k; i < m; i++) {
172 s += QR_(i, k) * Q(i, j);
173 }
174 s = -s / QR_(k, k);
175 for (i = k; i < m; i++) {
176 Q(i, j) += s * QR_(i, k);
177 }
178 }
179 }
180 }
181 return Q;
182 }
183
184 /** Least squares solution of A*x = b
185 @param b m-length array (vector).
186 @return x n-length array (vector) that minimizes the two norm of
187 Q*R*X-B. If B is non-conformant, or if QR.isFullRank() is false, the
188 routine returns a null (0-length) vector.
189 */
191 if (b.size() != (unsigned int)m) /* arrays must be conformant */
192 return DynamicVector<Real>();
193
194 if (!isFullRank()) /* matrix is rank deficient */
195 {
196 return DynamicVector<Real>();
197 }
198
200
201 // Compute Y = transpose(Q)*b
202 for (int k = 0; k < n; k++) {
203 Real s = 0.0;
204 for (int i = k; i < m; i++) {
205 s += QR_(i, k) * x(i);
206 }
207 s = -s / QR_(k, k);
208 for (int i = k; i < m; i++) {
209 x(i) += s * QR_(i, k);
210 }
211 }
212 // Solve R*X = Y;
213 for (int k = n - 1; k >= 0; k--) {
214 x(k) /= Rdiag(k);
215 for (int i = 0; i < k; i++) {
216 x(i) -= x(k) * QR_(i, k);
217 }
218 }
219
220 /* return n x nx portion of X */
222 for (int i = 0; i < n; i++)
223 x_(i) = x(i);
224
225 return x_;
226 }
227
228 /** Least squares solution of A*X = B
229 @param B m x k Array (must conform).
230 @return X n x k Array that minimizes the two norm of Q*R*X-B. If
231 B is non-conformant, or if QR.isFullRank() is false,
232 the routine returns a null (0x0) array.
233 */
235 if (B.getNRow() != m) /* arrays must be conformant */
236 return DynamicRectMatrix<Real>(0, 0);
237
238 if (!isFullRank()) /* matrix is rank deficient */
239 {
240 return DynamicRectMatrix<Real>(0, 0);
241 }
242
243 int nx = B.getNCol();
245 int i = 0, j = 0, k = 0;
246
247 // Compute Y = transpose(Q)*B
248 for (k = 0; k < n; k++) {
249 for (j = 0; j < nx; j++) {
250 Real s = 0.0;
251 for (i = k; i < m; i++) {
252 s += QR_(i, k) * X(i, j);
253 }
254 s = -s / QR_(k, k);
255 for (i = k; i < m; i++) {
256 X(i, j) += s * QR_(i, k);
257 }
258 }
259 }
260 // Solve R*X = Y;
261 for (k = n - 1; k >= 0; k--) {
262 for (j = 0; j < nx; j++) {
263 X(k, j) /= Rdiag(k);
264 }
265 for (i = 0; i < k; i++) {
266 for (j = 0; j < nx; j++) {
267 X(i, j) -= X(k, j) * QR_(i, k);
268 }
269 }
270 }
271
272 /* return n x nx portion of X */
273 DynamicRectMatrix<Real> X_(n, nx);
274 for (i = 0; i < n; i++)
275 for (j = 0; j < nx; j++)
276 X_(i, j) = X(i, j);
277
278 return X_;
279 }
280 };
281
282} // namespace JAMA
283// namespace JAMA
284
285#endif
286// JAMA_QR__H
DynamicRectMatrix< Real > solve(const DynamicRectMatrix< Real > &B) const
Least squares solution of A*X = B.
Definition QR.hpp:234
int isFullRank() const
Flag to denote the matrix is of full rank.
Definition QR.hpp:106
DynamicVector< Real > solve(const DynamicVector< Real > &b) const
Least squares solution of A*x = b.
Definition QR.hpp:190
DynamicRectMatrix< Real > getQ() const
Generate and return the (economy-sized) orthogonal factor.
Definition QR.hpp:159
DynamicRectMatrix< Real > getHouseholder(void) const
Retreive the Householder vectors from QR factorization.
Definition QR.hpp:118
DynamicRectMatrix< Real > getR() const
Return the upper triangular factor, R, of the QR factorization.
Definition QR.hpp:139
rectangular matrix class
Dynamically-sized vector class.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.