SingularValueDecomposition.java

  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.  * This is not the original file distributed by the Apache Software Foundation
  19.  * It has been modified by the Hipparchus project
  20.  */
  21. package org.hipparchus.linear;

  22. import org.hipparchus.exception.LocalizedCoreFormats;
  23. import org.hipparchus.exception.MathIllegalArgumentException;
  24. import org.hipparchus.util.FastMath;
  25. import org.hipparchus.util.Precision;

  26. /**
  27.  * Calculates the compact Singular Value Decomposition of a matrix.
  28.  * <p>
  29.  * The Singular Value Decomposition of matrix A is a set of three matrices: U,
  30.  * &Sigma; and V such that A = U &times; &Sigma; &times; V<sup>T</sup>. Let A be
  31.  * a m &times; n matrix, then U is a m &times; p orthogonal matrix, &Sigma; is a
  32.  * p &times; p diagonal matrix with positive or null elements, V is a p &times;
  33.  * n orthogonal matrix (hence V<sup>T</sup> is also orthogonal) where
  34.  * p=min(m,n).
  35.  * </p>
  36.  * <p>This class is similar to the class with similar name from the
  37.  * <a href="http://math.nist.gov/javanumerics/jama/">JAMA</a> library, with the
  38.  * following changes:</p>
  39.  * <ul>
  40.  *   <li>the {@code norm2} method which has been renamed as {@link #getNorm()
  41.  *   getNorm},</li>
  42.  *   <li>the {@code cond} method which has been renamed as {@link
  43.  *   #getConditionNumber() getConditionNumber},</li>
  44.  *   <li>the {@code rank} method which has been renamed as {@link #getRank()
  45.  *   getRank},</li>
  46.  *   <li>a {@link #getUT() getUT} method has been added,</li>
  47.  *   <li>a {@link #getVT() getVT} method has been added,</li>
  48.  *   <li>a {@link #getSolver() getSolver} method has been added,</li>
  49.  *   <li>a {@link #getCovariance(double) getCovariance} method has been added.</li>
  50.  * </ul>
  51.  * @see <a href="http://mathworld.wolfram.com/SingularValueDecomposition.html">MathWorld</a>
  52.  * @see <a href="http://en.wikipedia.org/wiki/Singular_value_decomposition">Wikipedia</a>
  53.  */
  54. public class SingularValueDecomposition {
  55.     /** Relative threshold for small singular values. */
  56.     private static final double EPS = 0x1.0p-52;
  57.     /** Absolute threshold for small singular values. */
  58.     private static final double TINY = 0x1.0p-966;
  59.     /** Computed singular values. */
  60.     private final double[] singularValues;
  61.     /** max(row dimension, column dimension). */
  62.     private final int m;
  63.     /** min(row dimension, column dimension). */
  64.     private final int n;
  65.     /** Cached value of U matrix. */
  66.     private final RealMatrix cachedU;
  67.     /** Cached value of transposed U matrix. */
  68.     private RealMatrix cachedUt;
  69.     /** Cached value of S (diagonal) matrix. */
  70.     private RealMatrix cachedS;
  71.     /** Cached value of V matrix. */
  72.     private final RealMatrix cachedV;
  73.     /** Cached value of transposed V matrix. */
  74.     private RealMatrix cachedVt;
  75.     /**
  76.      * Tolerance value for small singular values, calculated once we have
  77.      * populated "singularValues".
  78.      **/
  79.     private final double tol;

  80.     /**
  81.      * Calculates the compact Singular Value Decomposition of the given matrix.
  82.      *
  83.      * @param matrix Matrix to decompose.
  84.      */
  85.     public SingularValueDecomposition(final RealMatrix matrix) {
  86.         final double[][] A;

  87.          // "m" is always the largest dimension.
  88.         final boolean transposed;
  89.         if (matrix.getRowDimension() < matrix.getColumnDimension()) {
  90.             transposed = true;
  91.             A = matrix.transpose().getData();
  92.             m = matrix.getColumnDimension();
  93.             n = matrix.getRowDimension();
  94.         } else {
  95.             transposed = false;
  96.             A = matrix.getData();
  97.             m = matrix.getRowDimension();
  98.             n = matrix.getColumnDimension();
  99.         }

  100.         singularValues = new double[n];
  101.         final double[][] U = new double[m][n];
  102.         final double[][] V = new double[n][n];
  103.         final double[] e = new double[n];
  104.         final double[] work = new double[m];
  105.         // Reduce A to bidiagonal form, storing the diagonal elements
  106.         // in s and the super-diagonal elements in e.
  107.         final int nct = FastMath.min(m - 1, n);
  108.         final int nrt = FastMath.max(0, n - 2);
  109.         for (int k = 0; k < FastMath.max(nct, nrt); k++) {
  110.             if (k < nct) {
  111.                 // Compute the transformation for the k-th column and
  112.                 // place the k-th diagonal in s[k].
  113.                 // Compute 2-norm of k-th column without under/overflow.
  114.                 singularValues[k] = 0;
  115.                 for (int i = k; i < m; i++) {
  116.                     singularValues[k] = FastMath.hypot(singularValues[k], A[i][k]);
  117.                 }
  118.                 if (singularValues[k] != 0) {
  119.                     if (A[k][k] < 0) {
  120.                         singularValues[k] = -singularValues[k];
  121.                     }
  122.                     for (int i = k; i < m; i++) {
  123.                         A[i][k] /= singularValues[k];
  124.                     }
  125.                     A[k][k] += 1;
  126.                 }
  127.                 singularValues[k] = -singularValues[k];
  128.             }
  129.             for (int j = k + 1; j < n; j++) {
  130.                 if (k < nct &&
  131.                     singularValues[k] != 0) {
  132.                     // Apply the transformation.
  133.                     double t = 0;
  134.                     for (int i = k; i < m; i++) {
  135.                         t += A[i][k] * A[i][j];
  136.                     }
  137.                     t = -t / A[k][k];
  138.                     for (int i = k; i < m; i++) {
  139.                         A[i][j] += t * A[i][k];
  140.                     }
  141.                 }
  142.                 // Place the k-th row of A into e for the
  143.                 // subsequent calculation of the row transformation.
  144.                 e[j] = A[k][j];
  145.             }
  146.             if (k < nct) {
  147.                 // Place the transformation in U for subsequent back
  148.                 // multiplication.
  149.                 for (int i = k; i < m; i++) {
  150.                     U[i][k] = A[i][k];
  151.                 }
  152.             }
  153.             if (k < nrt) {
  154.                 // Compute the k-th row transformation and place the
  155.                 // k-th super-diagonal in e[k].
  156.                 // Compute 2-norm without under/overflow.
  157.                 e[k] = 0;
  158.                 for (int i = k + 1; i < n; i++) {
  159.                     e[k] = FastMath.hypot(e[k], e[i]);
  160.                 }
  161.                 if (e[k] != 0) {
  162.                     if (e[k + 1] < 0) {
  163.                         e[k] = -e[k];
  164.                     }
  165.                     for (int i = k + 1; i < n; i++) {
  166.                         e[i] /= e[k];
  167.                     }
  168.                     e[k + 1] += 1;
  169.                 }
  170.                 e[k] = -e[k];
  171.                 if (k + 1 < m &&
  172.                     e[k] != 0) {
  173.                     // Apply the transformation.
  174.                     for (int i = k + 1; i < m; i++) {
  175.                         work[i] = 0;
  176.                     }
  177.                     for (int j = k + 1; j < n; j++) {
  178.                         for (int i = k + 1; i < m; i++) {
  179.                             work[i] += e[j] * A[i][j];
  180.                         }
  181.                     }
  182.                     for (int j = k + 1; j < n; j++) {
  183.                         final double t = -e[j] / e[k + 1];
  184.                         for (int i = k + 1; i < m; i++) {
  185.                             A[i][j] += t * work[i];
  186.                         }
  187.                     }
  188.                 }

  189.                 // Place the transformation in V for subsequent
  190.                 // back multiplication.
  191.                 for (int i = k + 1; i < n; i++) {
  192.                     V[i][k] = e[i];
  193.                 }
  194.             }
  195.         }
  196.         // Set up the final bidiagonal matrix or order p.
  197.         int p = n;
  198.         if (nct < n) {
  199.             singularValues[nct] = A[nct][nct];
  200.         }
  201.         if (m < p) {
  202.             singularValues[p - 1] = 0;
  203.         }
  204.         if (nrt + 1 < p) {
  205.             e[nrt] = A[nrt][p - 1];
  206.         }
  207.         e[p - 1] = 0;

  208.         // Generate U.
  209.         for (int j = nct; j < n; j++) {
  210.             for (int i = 0; i < m; i++) {
  211.                 U[i][j] = 0;
  212.             }
  213.             U[j][j] = 1;
  214.         }
  215.         for (int k = nct - 1; k >= 0; k--) {
  216.             if (singularValues[k] != 0) {
  217.                 for (int j = k + 1; j < n; j++) {
  218.                     double t = 0;
  219.                     for (int i = k; i < m; i++) {
  220.                         t += U[i][k] * U[i][j];
  221.                     }
  222.                     t = -t / U[k][k];
  223.                     for (int i = k; i < m; i++) {
  224.                         U[i][j] += t * U[i][k];
  225.                     }
  226.                 }
  227.                 for (int i = k; i < m; i++) {
  228.                     U[i][k] = -U[i][k];
  229.                 }
  230.                 U[k][k] = 1 + U[k][k];
  231.                 for (int i = 0; i < k - 1; i++) {
  232.                     U[i][k] = 0;
  233.                 }
  234.             } else {
  235.                 for (int i = 0; i < m; i++) {
  236.                     U[i][k] = 0;
  237.                 }
  238.                 U[k][k] = 1;
  239.             }
  240.         }

  241.         // Generate V.
  242.         for (int k = n - 1; k >= 0; k--) {
  243.             if (k < nrt &&
  244.                 e[k] != 0) {
  245.                 for (int j = k + 1; j < n; j++) {
  246.                     double t = 0;
  247.                     for (int i = k + 1; i < n; i++) {
  248.                         t += V[i][k] * V[i][j];
  249.                     }
  250.                     t = -t / V[k + 1][k];
  251.                     for (int i = k + 1; i < n; i++) {
  252.                         V[i][j] += t * V[i][k];
  253.                     }
  254.                 }
  255.             }
  256.             for (int i = 0; i < n; i++) {
  257.                 V[i][k] = 0;
  258.             }
  259.             V[k][k] = 1;
  260.         }

  261.         // Main iteration loop for the singular values.
  262.         final int pp = p - 1;
  263.         while (p > 0) {
  264.             int k;
  265.             int kase;
  266.             // Here is where a test for too many iterations would go.
  267.             // This section of the program inspects for
  268.             // negligible elements in the s and e arrays.  On
  269.             // completion the variables kase and k are set as follows.
  270.             // kase = 1     if s(p) and e[k-1] are negligible and k<p
  271.             // kase = 2     if s(k) is negligible and k<p
  272.             // kase = 3     if e[k-1] is negligible, k<p, and
  273.             //              s(k), ..., s(p) are not negligible (qr step).
  274.             // kase = 4     if e(p-1) is negligible (convergence).
  275.             for (k = p - 2; k >= 0; k--) {
  276.                 final double threshold
  277.                     = TINY + EPS * (FastMath.abs(singularValues[k]) +
  278.                                     FastMath.abs(singularValues[k + 1]));

  279.                 // the following condition is written this way in order
  280.                 // to break out of the loop when NaN occurs, writing it
  281.                 // as "if (FastMath.abs(e[k]) <= threshold)" would loop
  282.                 // indefinitely in case of NaNs because comparison on NaNs
  283.                 // always return false, regardless of what is checked
  284.                 // see issue MATH-947
  285.                 if (!(FastMath.abs(e[k]) > threshold)) { // NOPMD - as explained above, the way this test is written is correct
  286.                     e[k] = 0;
  287.                     break;
  288.                 }

  289.             }

  290.             if (k == p - 2) {
  291.                 kase = 4;
  292.             } else {
  293.                 int ks;
  294.                 for (ks = p - 1; ks >= k; ks--) {
  295.                     if (ks == k) {
  296.                         break;
  297.                     }
  298.                     final double t = (ks != p ? FastMath.abs(e[ks]) : 0) +
  299.                         (ks != k + 1 ? FastMath.abs(e[ks - 1]) : 0);
  300.                     if (FastMath.abs(singularValues[ks]) <= TINY + EPS * t) {
  301.                         singularValues[ks] = 0;
  302.                         break;
  303.                     }
  304.                 }
  305.                 if (ks == k) {
  306.                     kase = 3;
  307.                 } else if (ks == p - 1) {
  308.                     kase = 1;
  309.                 } else {
  310.                     kase = 2;
  311.                     k = ks;
  312.                 }
  313.             }
  314.             k++;
  315.             // Perform the task indicated by kase.
  316.             switch (kase) { // NOPMD - breaking this complex algorithm into functions just to keep PMD happy would be artificial
  317.                 // Deflate negligible s(p).
  318.                 case 1: {
  319.                     double f = e[p - 2];
  320.                     e[p - 2] = 0;
  321.                     for (int j = p - 2; j >= k; j--) {
  322.                         double t = FastMath.hypot(singularValues[j], f);
  323.                         final double cs = singularValues[j] / t;
  324.                         final double sn = f / t;
  325.                         singularValues[j] = t;
  326.                         if (j != k) {
  327.                             f = -sn * e[j - 1];
  328.                             e[j - 1] = cs * e[j - 1];
  329.                         }

  330.                         for (int i = 0; i < n; i++) {
  331.                             t = cs * V[i][j] + sn * V[i][p - 1];
  332.                             V[i][p - 1] = -sn * V[i][j] + cs * V[i][p - 1];
  333.                             V[i][j] = t;
  334.                         }
  335.                     }
  336.                 }
  337.                 break;
  338.                 // Split at negligible s(k).
  339.                 case 2: {
  340.                     double f = e[k - 1];
  341.                     e[k - 1] = 0;
  342.                     for (int j = k; j < p; j++) {
  343.                         double t = FastMath.hypot(singularValues[j], f);
  344.                         final double cs = singularValues[j] / t;
  345.                         final double sn = f / t;
  346.                         singularValues[j] = t;
  347.                         f = -sn * e[j];
  348.                         e[j] = cs * e[j];

  349.                         for (int i = 0; i < m; i++) {
  350.                             t = cs * U[i][j] + sn * U[i][k - 1];
  351.                             U[i][k - 1] = -sn * U[i][j] + cs * U[i][k - 1];
  352.                             U[i][j] = t;
  353.                         }
  354.                     }
  355.                 }
  356.                 break;
  357.                 // Perform one qr step.
  358.                 case 3: {
  359.                     // Calculate the shift.
  360.                     final double maxPm1Pm2 = FastMath.max(FastMath.abs(singularValues[p - 1]),
  361.                                                           FastMath.abs(singularValues[p - 2]));
  362.                     final double scale = FastMath.max(FastMath.max(FastMath.max(maxPm1Pm2,
  363.                                                                                 FastMath.abs(e[p - 2])),
  364.                                                                    FastMath.abs(singularValues[k])),
  365.                                                       FastMath.abs(e[k]));
  366.                     final double sp = singularValues[p - 1] / scale;
  367.                     final double spm1 = singularValues[p - 2] / scale;
  368.                     final double epm1 = e[p - 2] / scale;
  369.                     final double sk = singularValues[k] / scale;
  370.                     final double ek = e[k] / scale;
  371.                     final double b = ((spm1 + sp) * (spm1 - sp) + epm1 * epm1) / 2.0;
  372.                     final double c = (sp * epm1) * (sp * epm1);
  373.                     double shift = 0;
  374.                     if (b != 0 ||
  375.                         c != 0) {
  376.                         shift = FastMath.sqrt(b * b + c);
  377.                         if (b < 0) {
  378.                             shift = -shift;
  379.                         }
  380.                         shift = c / (b + shift);
  381.                     }
  382.                     double f = (sk + sp) * (sk - sp) + shift;
  383.                     double g = sk * ek;
  384.                     // Chase zeros.
  385.                     for (int j = k; j < p - 1; j++) {
  386.                         double t = FastMath.hypot(f, g);
  387.                         double cs = f / t;
  388.                         double sn = g / t;
  389.                         if (j != k) {
  390.                             e[j - 1] = t;
  391.                         }
  392.                         f = cs * singularValues[j] + sn * e[j];
  393.                         e[j] = cs * e[j] - sn * singularValues[j];
  394.                         g = sn * singularValues[j + 1];
  395.                         singularValues[j + 1] = cs * singularValues[j + 1];

  396.                         for (int i = 0; i < n; i++) {
  397.                             t = cs * V[i][j] + sn * V[i][j + 1];
  398.                             V[i][j + 1] = -sn * V[i][j] + cs * V[i][j + 1];
  399.                             V[i][j] = t;
  400.                         }
  401.                         t = FastMath.hypot(f, g);
  402.                         cs = f / t;
  403.                         sn = g / t;
  404.                         singularValues[j] = t;
  405.                         f = cs * e[j] + sn * singularValues[j + 1];
  406.                         singularValues[j + 1] = -sn * e[j] + cs * singularValues[j + 1];
  407.                         g = sn * e[j + 1];
  408.                         e[j + 1] = cs * e[j + 1];
  409.                         if (j < m - 1) {
  410.                             for (int i = 0; i < m; i++) {
  411.                                 t = cs * U[i][j] + sn * U[i][j + 1];
  412.                                 U[i][j + 1] = -sn * U[i][j] + cs * U[i][j + 1];
  413.                                 U[i][j] = t;
  414.                             }
  415.                         }
  416.                     }
  417.                     e[p - 2] = f;
  418.                 }
  419.                 break;
  420.                 // Convergence.
  421.                 default: {
  422.                     // Make the singular values positive.
  423.                     if (singularValues[k] <= 0) {
  424.                         singularValues[k] = singularValues[k] < 0 ? -singularValues[k] : 0;

  425.                         for (int i = 0; i <= pp; i++) {
  426.                             V[i][k] = -V[i][k];
  427.                         }
  428.                     }
  429.                     // Order the singular values.
  430.                     while (k < pp) {
  431.                         if (singularValues[k] >= singularValues[k + 1]) {
  432.                             break;
  433.                         }
  434.                         double t = singularValues[k];
  435.                         singularValues[k] = singularValues[k + 1];
  436.                         singularValues[k + 1] = t;
  437.                         if (k < n - 1) {
  438.                             for (int i = 0; i < n; i++) {
  439.                                 t = V[i][k + 1];
  440.                                 V[i][k + 1] = V[i][k];
  441.                                 V[i][k] = t;
  442.                             }
  443.                         }
  444.                         if (k < m - 1) {
  445.                             for (int i = 0; i < m; i++) {
  446.                                 t = U[i][k + 1];
  447.                                 U[i][k + 1] = U[i][k];
  448.                                 U[i][k] = t;
  449.                             }
  450.                         }
  451.                         k++;
  452.                     }
  453.                     p--;
  454.                 }
  455.                 break;
  456.             }
  457.         }

  458.         // Set the small value tolerance used to calculate rank and pseudo-inverse
  459.         tol = FastMath.max(m * singularValues[0] * EPS,
  460.                            FastMath.sqrt(Precision.SAFE_MIN));

  461.         if (!transposed) {
  462.             cachedU = MatrixUtils.createRealMatrix(U);
  463.             cachedV = MatrixUtils.createRealMatrix(V);
  464.         } else {
  465.             cachedU = MatrixUtils.createRealMatrix(V);
  466.             cachedV = MatrixUtils.createRealMatrix(U);
  467.         }
  468.     }

  469.     /**
  470.      * Returns the matrix U of the decomposition.
  471.      * <p>U is an orthogonal matrix, i.e. its transpose is also its inverse.</p>
  472.      * @return the U matrix
  473.      * @see #getUT()
  474.      */
  475.     public RealMatrix getU() {
  476.         // return the cached matrix
  477.         return cachedU;

  478.     }

  479.     /**
  480.      * Returns the transpose of the matrix U of the decomposition.
  481.      * <p>U is an orthogonal matrix, i.e. its transpose is also its inverse.</p>
  482.      * @return the U matrix (or null if decomposed matrix is singular)
  483.      * @see #getU()
  484.      */
  485.     public RealMatrix getUT() {
  486.         if (cachedUt == null) {
  487.             cachedUt = getU().transpose();
  488.         }
  489.         // return the cached matrix
  490.         return cachedUt;
  491.     }

  492.     /**
  493.      * Returns the diagonal matrix &Sigma; of the decomposition.
  494.      * <p>&Sigma; is a diagonal matrix. The singular values are provided in
  495.      * non-increasing order, for compatibility with Jama.</p>
  496.      * @return the &Sigma; matrix
  497.      */
  498.     public RealMatrix getS() {
  499.         if (cachedS == null) {
  500.             // cache the matrix for subsequent calls
  501.             cachedS = MatrixUtils.createRealDiagonalMatrix(singularValues);
  502.         }
  503.         return cachedS;
  504.     }

  505.     /**
  506.      * Returns the diagonal elements of the matrix &Sigma; of the decomposition.
  507.      * <p>The singular values are provided in non-increasing order, for
  508.      * compatibility with Jama.</p>
  509.      * @return the diagonal elements of the &Sigma; matrix
  510.      */
  511.     public double[] getSingularValues() {
  512.         return singularValues.clone();
  513.     }

  514.     /**
  515.      * Returns the matrix V of the decomposition.
  516.      * <p>V is an orthogonal matrix, i.e. its transpose is also its inverse.</p>
  517.      * @return the V matrix (or null if decomposed matrix is singular)
  518.      * @see #getVT()
  519.      */
  520.     public RealMatrix getV() {
  521.         // return the cached matrix
  522.         return cachedV;
  523.     }

  524.     /**
  525.      * Returns the transpose of the matrix V of the decomposition.
  526.      * <p>V is an orthogonal matrix, i.e. its transpose is also its inverse.</p>
  527.      * @return the V matrix (or null if decomposed matrix is singular)
  528.      * @see #getV()
  529.      */
  530.     public RealMatrix getVT() {
  531.         if (cachedVt == null) {
  532.             cachedVt = getV().transpose();
  533.         }
  534.         // return the cached matrix
  535.         return cachedVt;
  536.     }

  537.     /**
  538.      * Returns the n &times; n covariance matrix.
  539.      * <p>The covariance matrix is V &times; J &times; V<sup>T</sup>
  540.      * where J is the diagonal matrix of the inverse of the squares of
  541.      * the singular values.</p>
  542.      * @param minSingularValue value below which singular values are ignored
  543.      * (a 0 or negative value implies all singular value will be used)
  544.      * @return covariance matrix
  545.      * @exception IllegalArgumentException if minSingularValue is larger than
  546.      * the largest singular value, meaning all singular values are ignored
  547.      */
  548.     public RealMatrix getCovariance(final double minSingularValue) {
  549.         // get the number of singular values to consider
  550.         final int p = singularValues.length;
  551.         int dimension = 0;
  552.         while (dimension < p &&
  553.                singularValues[dimension] >= minSingularValue) {
  554.             ++dimension;
  555.         }

  556.         if (dimension == 0) {
  557.             throw new MathIllegalArgumentException(LocalizedCoreFormats.TOO_LARGE_CUTOFF_SINGULAR_VALUE,
  558.                                                 minSingularValue, singularValues[0], true);
  559.         }

  560.         final double[][] data = new double[dimension][p];
  561.         getVT().walkInOptimizedOrder(new DefaultRealMatrixPreservingVisitor() {
  562.             /** {@inheritDoc} */
  563.             @Override
  564.             public void visit(final int row, final int column,
  565.                     final double value) {
  566.                 data[row][column] = value / singularValues[row];
  567.             }
  568.         }, 0, dimension - 1, 0, p - 1);

  569.         RealMatrix jv = new Array2DRowRealMatrix(data, false);
  570.         return jv.transposeMultiply(jv);
  571.     }

  572.     /**
  573.      * Returns the L<sub>2</sub> norm of the matrix.
  574.      * <p>The L<sub>2</sub> norm is max(|A &times; u|<sub>2</sub> /
  575.      * |u|<sub>2</sub>), where |.|<sub>2</sub> denotes the vectorial 2-norm
  576.      * (i.e. the traditional euclidian norm).</p>
  577.      * @return norm
  578.      */
  579.     public double getNorm() {
  580.         return singularValues[0];
  581.     }

  582.     /**
  583.      * Return the condition number of the matrix.
  584.      * @return condition number of the matrix
  585.      */
  586.     public double getConditionNumber() {
  587.         return singularValues[0] / singularValues[n - 1];
  588.     }

  589.     /**
  590.      * Computes the inverse of the condition number.
  591.      * In cases of rank deficiency, the {@link #getConditionNumber() condition
  592.      * number} will become undefined.
  593.      *
  594.      * @return the inverse of the condition number.
  595.      */
  596.     public double getInverseConditionNumber() {
  597.         return singularValues[n - 1] / singularValues[0];
  598.     }

  599.     /**
  600.      * Return the effective numerical matrix rank.
  601.      * <p>The effective numerical rank is the number of non-negligible
  602.      * singular values. The threshold used to identify non-negligible
  603.      * terms is max(m,n) &times; ulp(s<sub>1</sub>) where ulp(s<sub>1</sub>)
  604.      * is the least significant bit of the largest singular value.</p>
  605.      * @return effective numerical matrix rank
  606.      */
  607.     public int getRank() {
  608.         int r = 0;
  609.         for (int i = 0; i < singularValues.length; i++) {
  610.             if (singularValues[i] > tol) {
  611.                 r++;
  612.             }
  613.         }
  614.         return r;
  615.     }

  616.     /**
  617.      * Get a solver for finding the A &times; X = B solution in least square sense.
  618.      * @return a solver
  619.      */
  620.     public DecompositionSolver getSolver() {
  621.         return new Solver(singularValues, getUT(), getV(), getRank() == m, tol);
  622.     }

  623.     /** Specialized solver. */
  624.     private static class Solver implements DecompositionSolver {
  625.         /** Pseudo-inverse of the initial matrix. */
  626.         private final RealMatrix pseudoInverse;
  627.         /** Singularity indicator. */
  628.         private final boolean nonSingular;

  629.         /**
  630.          * Build a solver from decomposed matrix.
  631.          *
  632.          * @param singularValues Singular values.
  633.          * @param uT U<sup>T</sup> matrix of the decomposition.
  634.          * @param v V matrix of the decomposition.
  635.          * @param nonSingular Singularity indicator.
  636.          * @param tol tolerance for singular values
  637.          */
  638.         private Solver(final double[] singularValues, final RealMatrix uT,
  639.                        final RealMatrix v, final boolean nonSingular, final double tol) {
  640.             final double[][] suT = uT.getData();
  641.             for (int i = 0; i < singularValues.length; ++i) {
  642.                 final double a;
  643.                 if (singularValues[i] > tol) {
  644.                     a = 1 / singularValues[i];
  645.                 } else {
  646.                     a = 0;
  647.                 }
  648.                 final double[] suTi = suT[i];
  649.                 for (int j = 0; j < suTi.length; ++j) {
  650.                     suTi[j] *= a;
  651.                 }
  652.             }
  653.             pseudoInverse = v.multiply(new Array2DRowRealMatrix(suT, false));
  654.             this.nonSingular = nonSingular;
  655.         }

  656.         /**
  657.          * Solve the linear equation A &times; X = B in least square sense.
  658.          * <p>
  659.          * The m&times;n matrix A may not be square, the solution X is such that
  660.          * ||A &times; X - B|| is minimal.
  661.          * </p>
  662.          * @param b Right-hand side of the equation A &times; X = B
  663.          * @return a vector X that minimizes the two norm of A &times; X - B
  664.          * @throws org.hipparchus.exception.MathIllegalArgumentException
  665.          * if the matrices dimensions do not match.
  666.          */
  667.         @Override
  668.         public RealVector solve(final RealVector b) {
  669.             return pseudoInverse.operate(b);
  670.         }

  671.         /**
  672.          * Solve the linear equation A &times; X = B in least square sense.
  673.          * <p>
  674.          * The m&times;n matrix A may not be square, the solution X is such that
  675.          * ||A &times; X - B|| is minimal.
  676.          * </p>
  677.          *
  678.          * @param b Right-hand side of the equation A &times; X = B
  679.          * @return a matrix X that minimizes the two norm of A &times; X - B
  680.          * @throws org.hipparchus.exception.MathIllegalArgumentException
  681.          * if the matrices dimensions do not match.
  682.          */
  683.         @Override
  684.         public RealMatrix solve(final RealMatrix b) {
  685.             return pseudoInverse.multiply(b);
  686.         }

  687.         /**
  688.          * Check if the decomposed matrix is non-singular.
  689.          *
  690.          * @return {@code true} if the decomposed matrix is non-singular.
  691.          */
  692.         @Override
  693.         public boolean isNonSingular() {
  694.             return nonSingular;
  695.         }

  696.         /**
  697.          * Get the pseudo-inverse of the decomposed matrix.
  698.          *
  699.          * @return the inverse matrix.
  700.          */
  701.         @Override
  702.         public RealMatrix getInverse() {
  703.             return pseudoInverse;
  704.         }

  705.         /** {@inheritDoc} */
  706.         @Override
  707.         public int getRowDimension() {
  708.             return pseudoInverse.getColumnDimension();
  709.         }

  710.         /** {@inheritDoc} */
  711.         @Override
  712.         public int getColumnDimension() {
  713.             return pseudoInverse.getRowDimension();
  714.         }

  715.     }
  716. }