1 /*
2 * Licensed to the Apache Software Foundation (ASF) under one or more
3 * contributor license agreements. See the NOTICE file distributed with
4 * this work for additional information regarding copyright ownership.
5 * The ASF licenses this file to You under the Apache License, Version 2.0
6 * (the "License"); you may not use this file except in compliance with
7 * the License. You may obtain a copy of the License at
8 *
9 * https://www.apache.org/licenses/LICENSE-2.0
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 */
17
18 /*
19 * This is not the original file distributed by the Apache Software Foundation
20 * It has been modified by the Hipparchus project
21 */
22
23 package org.hipparchus.linear;
24
25 import org.hipparchus.exception.LocalizedCoreFormats;
26 import org.hipparchus.exception.MathIllegalArgumentException;
27 import org.hipparchus.util.FastMath;
28
29 /**
30 * Calculates the rectangular Cholesky decomposition of a matrix.
31 * <p>The rectangular Cholesky decomposition of a real symmetric positive
32 * semidefinite matrix A consists of a rectangular matrix B with the same
33 * number of rows such that: A is almost equal to BB<sup>T</sup>, depending
34 * on a user-defined tolerance. In a sense, this is the square root of A.</p>
35 * <p>The difference with respect to the regular {@link CholeskyDecomposition}
36 * is that rows/columns may be permuted (hence the rectangular shape instead
37 * of the traditional triangular shape) and there is a threshold to ignore
38 * small diagonal elements. This is used for example to generate {@link
39 * org.hipparchus.random.CorrelatedRandomVectorGenerator correlated
40 * random n-dimensions vectors} in a p-dimension subspace (p < n).
41 * In other words, it allows generating random vectors from a covariance
42 * matrix that is only positive semidefinite, and not positive definite.</p>
43 * <p>Rectangular Cholesky decomposition is <em>not</em> suited for solving
44 * linear systems, so it does not provide any {@link DecompositionSolver
45 * decomposition solver}.</p>
46 *
47 * @see <a href="http://mathworld.wolfram.com/CholeskyDecomposition.html">MathWorld</a>
48 * @see <a href="http://en.wikipedia.org/wiki/Cholesky_decomposition">Wikipedia</a>
49 */
50 public class RectangularCholeskyDecomposition {
51
52 /** Permutated Cholesky root of the symmetric positive semidefinite matrix. */
53 private final RealMatrix root;
54
55 /** Rank of the symmetric positive semidefinite matrix. */
56 private int rank;
57
58 /**
59 * Decompose a symmetric positive semidefinite matrix.
60 * <p>
61 * <b>Note:</b> this constructor follows the linpack method to detect dependent
62 * columns by proceeding with the Cholesky algorithm until a nonpositive diagonal
63 * element is encountered.
64 *
65 * @see <a href="http://eprints.ma.man.ac.uk/1193/01/covered/MIMS_ep2008_56.pdf">
66 * Analysis of the Cholesky Decomposition of a Semi-definite Matrix</a>
67 *
68 * @param matrix Symmetric positive semidefinite matrix.
69 * @exception MathIllegalArgumentException if the matrix is not
70 * positive semidefinite.
71 */
72 public RectangularCholeskyDecomposition(RealMatrix matrix)
73 throws MathIllegalArgumentException {
74 this(matrix, 0);
75 }
76
77 /**
78 * Decompose a symmetric positive semidefinite matrix.
79 *
80 * @param matrix Symmetric positive semidefinite matrix.
81 * @param small Diagonal elements threshold under which columns are
82 * considered to be dependent on previous ones and are discarded.
83 * @exception MathIllegalArgumentException if the matrix is not
84 * positive semidefinite.
85 */
86 public RectangularCholeskyDecomposition(RealMatrix matrix, double small)
87 throws MathIllegalArgumentException {
88
89 final int order = matrix.getRowDimension();
90 final double[][] c = matrix.getData();
91 final double[][] b = new double[order][order];
92
93 int[] index = new int[order];
94 for (int i = 0; i < order; ++i) {
95 index[i] = i;
96 }
97
98 int r = 0;
99 for (boolean loop = true; loop;) {
100
101 // find maximal diagonal element
102 int swapR = r;
103 for (int i = r + 1; i < order; ++i) {
104 int ii = index[i];
105 int isr = index[swapR];
106 if (c[ii][ii] > c[isr][isr]) {
107 swapR = i;
108 }
109 }
110
111
112 // swap elements
113 if (swapR != r) {
114 final int tmpIndex = index[r];
115 index[r] = index[swapR];
116 index[swapR] = tmpIndex;
117 final double[] tmpRow = b[r];
118 b[r] = b[swapR];
119 b[swapR] = tmpRow;
120 }
121
122 // check diagonal element
123 int ir = index[r];
124 if (c[ir][ir] <= small) {
125
126 if (r == 0) {
127 throw new MathIllegalArgumentException(LocalizedCoreFormats.NOT_POSITIVE_DEFINITE_MATRIX);
128 }
129
130 // check remaining diagonal elements
131 for (int i = r; i < order; ++i) {
132 if (c[index[i]][index[i]] < -small) {
133 // there is at least one sufficiently negative diagonal element,
134 // the symmetric positive semidefinite matrix is wrong
135 throw new MathIllegalArgumentException(LocalizedCoreFormats.NOT_POSITIVE_DEFINITE_MATRIX);
136 }
137 }
138
139 // all remaining diagonal elements are close to zero, we consider we have
140 // found the rank of the symmetric positive semidefinite matrix
141 loop = false;
142
143 } else {
144
145 // transform the matrix
146 final double sqrt = FastMath.sqrt(c[ir][ir]);
147 b[r][r] = sqrt;
148 final double inverse = 1 / sqrt;
149 final double inverse2 = 1 / c[ir][ir];
150 for (int i = r + 1; i < order; ++i) {
151 final int ii = index[i];
152 final double e = inverse * c[ii][ir];
153 b[i][r] = e;
154 c[ii][ii] -= c[ii][ir] * c[ii][ir] * inverse2;
155 for (int j = r + 1; j < i; ++j) {
156 final int ij = index[j];
157 final double f = c[ii][ij] - e * b[j][r];
158 c[ii][ij] = f;
159 c[ij][ii] = f;
160 }
161 }
162
163 // prepare next iteration
164 loop = ++r < order;
165 }
166 }
167
168 // build the root matrix
169 rank = r;
170 root = MatrixUtils.createRealMatrix(order, r);
171 for (int i = 0; i < order; ++i) {
172 for (int j = 0; j < r; ++j) {
173 root.setEntry(index[i], j, b[i][j]);
174 }
175 }
176
177 }
178
179 /** Get the root of the covariance matrix.
180 * The root is the rectangular matrix <code>B</code> such that
181 * the covariance matrix is equal to <code>B.B<sup>T</sup></code>
182 * @return root of the square matrix
183 * @see #getRank()
184 */
185 public RealMatrix getRootMatrix() {
186 return root;
187 }
188
189 /** Get the rank of the symmetric positive semidefinite matrix.
190 * The r is the number of independent rows in the symmetric positive semidefinite
191 * matrix, it is also the number of columns of the rectangular
192 * matrix of the decomposition.
193 * @return r of the square matrix.
194 * @see #getRootMatrix()
195 */
196 public int getRank() {
197 return rank;
198 }
199
200 }