summaryrefslogtreecommitdiff
path: root/src/de/lmu/ifi/dbs/elki/math/linearalgebra/QRDecomposition.java
blob: 5b52d837d0d235e541e6913c649549a69addd9c4 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
package de.lmu.ifi.dbs.elki.math.linearalgebra;

/*
 This file is part of ELKI:
 Environment for Developing KDD-Applications Supported by Index-Structures

 Copyright (C) 2012
 Ludwig-Maximilians-Universität München
 Lehr- und Forschungseinheit für Datenbanksysteme
 ELKI Development Team

 This program is free software: you can redistribute it and/or modify
 it under the terms of the GNU Affero General Public License as published by
 the Free Software Foundation, either version 3 of the License, or
 (at your option) any later version.

 This program is distributed in the hope that it will be useful,
 but WITHOUT ANY WARRANTY; without even the implied warranty of
 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 GNU Affero General Public License for more details.

 You should have received a copy of the GNU Affero General Public License
 along with this program.  If not, see <http://www.gnu.org/licenses/>.
 */

import de.lmu.ifi.dbs.elki.math.MathUtil;

/**
 * QR Decomposition.
 * <P>
 * For an m-by-n matrix A with m >= n, the QR decomposition is an m-by-n
 * orthogonal matrix Q and an n-by-n upper triangular matrix R so that A = Q*R.
 * <P>
 * The QR decompostion always exists, even if the matrix does not have full
 * rank, so the constructor will never fail. The primary use of the QR
 * decomposition is in the least squares solution of nonsquare systems of
 * simultaneous linear equations. This will fail if isFullRank() returns false.
 * 
 * @apiviz.uses Matrix - - transforms
 */
public class QRDecomposition implements java.io.Serializable {
  /**
   * Serial version
   */
  private static final long serialVersionUID = 1L;

  /**
   * Array for internal storage of decomposition.
   * 
   * @serial internal array storage.
   */
  private double[][] QR;

  /**
   * Row and column dimensions.
   * 
   * @serial column dimension.
   * @serial row dimension.
   */
  private int m, n;

  /**
   * Array for internal storage of diagonal of R.
   * 
   * @serial diagonal of R.
   */
  private double[] Rdiag;

  /*
   * ------------------------ Constructor ------------------------
   */

  /**
   * QR Decomposition, computed by Householder reflections.
   * 
   * @param A Rectangular matrix
   */
  public QRDecomposition(Matrix A) {
    this(A.getArrayRef(), A.getRowDimensionality(), A.getColumnDimensionality());
  }

  /**
   * QR Decomposition, computed by Householder reflections.
   * 
   * @param A Rectangular matrix
   * @param m row dimensionality
   * @param n column dimensionality
   */
  public QRDecomposition(double[][] A, int m, int n) {
    this.QR = new Matrix(A).getArrayCopy();
    this.m = QR.length;
    this.n = QR[0].length;
    Rdiag = new double[n];

    // Main loop.
    for(int k = 0; k < n; k++) {
      // Compute 2-norm of k-th column without under/overflow.
      double nrm = 0;
      for(int i = k; i < m; i++) {
        nrm = MathUtil.fastHypot(nrm, QR[i][k]);
      }

      if(nrm != 0.0) {
        // Form k-th Householder vector.
        if(QR[k][k] < 0) {
          nrm = -nrm;
        }
        for(int i = k; i < m; i++) {
          QR[i][k] /= nrm;
        }
        QR[k][k] += 1.0;

        // Apply transformation to remaining columns.
        for(int j = k + 1; j < n; j++) {
          double s = 0.0;
          for(int i = k; i < m; i++) {
            s += QR[i][k] * QR[i][j];
          }
          s = -s / QR[k][k];
          for(int i = k; i < m; i++) {
            QR[i][j] += s * QR[i][k];
          }
        }
      }
      Rdiag[k] = -nrm;
    }
  }

  /*
   * ------------------------ Public Methods ------------------------
   */

  /**
   * Is the matrix full rank?
   * 
   * @return true if R, and hence A, has full rank.
   */
  public boolean isFullRank() {
    for(int j = 0; j < n; j++) {
      if(Rdiag[j] == 0) {
        return false;
      }
    }
    return true;
  }

  /**
   * Return the Householder vectors
   * 
   * @return Lower trapezoidal matrix whose columns define the reflections
   */
  public Matrix getH() {
    Matrix X = new Matrix(m, n);
    double[][] H = X.getArrayRef();
    for(int i = 0; i < m; i++) {
      for(int j = 0; j < n; j++) {
        if(i >= j) {
          H[i][j] = QR[i][j];
        }
        else {
          H[i][j] = 0.0;
        }
      }
    }
    return X;
  }

  /**
   * Return the upper triangular factor
   * 
   * @return R
   */
  public Matrix getR() {
    Matrix X = new Matrix(n, n);
    double[][] R = X.getArrayRef();
    for(int i = 0; i < n; i++) {
      for(int j = 0; j < n; j++) {
        if(i < j) {
          R[i][j] = QR[i][j];
        }
        else if(i == j) {
          R[i][j] = Rdiag[i];
        }
        else {
          R[i][j] = 0.0;
        }
      }
    }
    return X;
  }

  /**
   * Generate and return the (economy-sized) orthogonal factor
   * 
   * @return Q
   */
  public Matrix getQ() {
    Matrix X = new Matrix(m, n);
    double[][] Q = X.getArrayRef();
    for(int k = n - 1; k >= 0; k--) {
      for(int i = 0; i < m; i++) {
        Q[i][k] = 0.0;
      }
      Q[k][k] = 1.0;
      for(int j = k; j < n; j++) {
        if(QR[k][k] != 0) {
          double s = 0.0;
          for(int i = k; i < m; i++) {
            s += QR[i][k] * Q[i][j];
          }
          s = -s / QR[k][k];
          for(int i = k; i < m; i++) {
            Q[i][j] += s * QR[i][k];
          }
        }
      }
    }
    return X;
  }

  /**
   * Least squares solution of A*X = B
   * 
   * @param B A Matrix with as many rows as A and any number of columns.
   * @return X that minimizes the two norm of Q*R*X-B.
   * @exception IllegalArgumentException Matrix row dimensions must agree.
   * @exception RuntimeException Matrix is rank deficient.
   */
  public Matrix solve(Matrix B) {
    if(B.getRowDimensionality() != m) {
      throw new IllegalArgumentException("Matrix row dimensions must agree.");
    }
    if(!this.isFullRank()) {
      throw new RuntimeException("Matrix is rank deficient.");
    }

    // Copy right hand side
    int nx = B.getColumnDimensionality();
    Matrix X = B.copy();

    solveInplace(X.getArrayRef(), nx);
    return X.getMatrix(0, n - 1, 0, nx - 1);
  }

  /**
   * Least squares solution of A*X = B
   * 
   * @param B A Matrix with as many rows as A and any number of columns.
   * @return X that minimizes the two norm of Q*R*X-B.
   * @exception IllegalArgumentException Matrix row dimensions must agree.
   * @exception RuntimeException Matrix is rank deficient.
   */
  public double[][] solve(double[][] B) {
    int rows = B.length;
    int cols = B[0].length;
    if(rows != m) {
      throw new IllegalArgumentException("Matrix row dimensions must agree.");
    }
    if(!this.isFullRank()) {
      throw new RuntimeException("Matrix is rank deficient.");
    }

    // Copy right hand side
    Matrix X = new Matrix(B).copy();

    solveInplace(X.getArrayRef(), cols);
    return X.getMatrix(0, n - 1, 0, cols - 1).getArrayRef();
  }

  private void solveInplace(double[][] X, int nx) {
    // Compute Y = transpose(Q)*B
    for(int k = 0; k < n; k++) {
      for(int j = 0; j < nx; j++) {
        double s = 0.0;
        for(int i = k; i < m; i++) {
          s += QR[i][k] * X[i][j];
        }
        s = -s / QR[k][k];
        for(int i = k; i < m; i++) {
          X[i][j] += s * QR[i][k];
        }
      }
    }
    // Solve R*X = Y;
    for(int k = n - 1; k >= 0; k--) {
      for(int j = 0; j < nx; j++) {
        X[k][j] /= Rdiag[k];
      }
      for(int i = 0; i < k; i++) {
        for(int j = 0; j < nx; j++) {
          X[i][j] -= X[k][j] * QR[i][k];
        }
      }
    }
  }
}