JacobiPreconditioner.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.analysis.function.Sqrt;
  23. import org.hipparchus.exception.LocalizedCoreFormats;
  24. import org.hipparchus.exception.MathIllegalArgumentException;
  25. import org.hipparchus.util.MathArrays;

  26. /**
  27.  * This class implements the standard Jacobi (diagonal) preconditioner. For a
  28.  * matrix A<sub>ij</sub>, this preconditioner is
  29.  * M = diag(1 / A<sub>11</sub>, 1 / A<sub>22</sub>, &hellip;).
  30.  */
  31. public class JacobiPreconditioner implements RealLinearOperator {

  32.     /** The diagonal coefficients of the preconditioner. */
  33.     private final ArrayRealVector diag;

  34.     /**
  35.      * Creates a new instance of this class.
  36.      *
  37.      * @param diag the diagonal coefficients of the linear operator to be
  38.      * preconditioned
  39.      * @param deep {@code true} if a deep copy of the above array should be
  40.      * performed
  41.      */
  42.     public JacobiPreconditioner(final double[] diag, final boolean deep) {
  43.         this.diag = new ArrayRealVector(diag, deep);
  44.     }

  45.     /**
  46.      * Creates a new instance of this class. This method extracts the diagonal
  47.      * coefficients of the specified linear operator. If {@code a} does not
  48.      * extend {@link AbstractRealMatrix}, then the coefficients of the
  49.      * underlying matrix are not accessible, coefficient extraction is made by
  50.      * matrix-vector products with the basis vectors (and might therefore take
  51.      * some time). With matrices, direct entry access is carried out.
  52.      *
  53.      * @param a the linear operator for which the preconditioner should be built
  54.      * @return the diagonal preconditioner made of the inverse of the diagonal
  55.      * coefficients of the specified linear operator
  56.      * @throws MathIllegalArgumentException if {@code a} is not square
  57.      */
  58.     public static JacobiPreconditioner create(final RealLinearOperator a)
  59.         throws MathIllegalArgumentException {
  60.         final int n = a.getColumnDimension();
  61.         if (a.getRowDimension() != n) {
  62.             throw new MathIllegalArgumentException(LocalizedCoreFormats.NON_SQUARE_OPERATOR,
  63.                                                    a.getRowDimension(), n);
  64.         }
  65.         final double[] diag = new double[n];
  66.         if (a instanceof AbstractRealMatrix) {
  67.             final AbstractRealMatrix m = (AbstractRealMatrix) a;
  68.             for (int i = 0; i < n; i++) {
  69.                 diag[i] = m.getEntry(i, i);
  70.             }
  71.         } else {
  72.             final ArrayRealVector x = new ArrayRealVector(n);
  73.             for (int i = 0; i < n; i++) {
  74.                 x.set(0.);
  75.                 x.setEntry(i, 1.);
  76.                 diag[i] = a.operate(x).getEntry(i);
  77.             }
  78.         }
  79.         return new JacobiPreconditioner(diag, false);
  80.     }

  81.     /** {@inheritDoc} */
  82.     @Override
  83.     public int getColumnDimension() {
  84.         return diag.getDimension();
  85.     }

  86.     /** {@inheritDoc} */
  87.     @Override
  88.     public int getRowDimension() {
  89.         return diag.getDimension();
  90.     }

  91.     /** {@inheritDoc} */
  92.     @Override
  93.     public RealVector operate(final RealVector x) {
  94.         // Dimension check is carried out by ebeDivide
  95.         return new ArrayRealVector(MathArrays.ebeDivide(x.toArray(),
  96.                                                         diag.toArray()),
  97.                                    false);
  98.     }

  99.     /**
  100.      * Returns the square root of {@code this} diagonal operator. More
  101.      * precisely, this method returns
  102.      * P = diag(1 / &radic;A<sub>11</sub>, 1 / &radic;A<sub>22</sub>, &hellip;).
  103.      *
  104.      * @return the square root of {@code this} preconditioner
  105.      */
  106.     public RealLinearOperator sqrt() {
  107.         final RealVector sqrtDiag = diag.map(new Sqrt());
  108.         return new RealLinearOperator() {
  109.             /** {@inheritDoc} */
  110.             @Override
  111.             public RealVector operate(final RealVector x) {
  112.                 return new ArrayRealVector(MathArrays.ebeDivide(x.toArray(),
  113.                                                                 sqrtDiag.toArray()),
  114.                                            false);
  115.             }

  116.             /** {@inheritDoc} */
  117.             @Override
  118.             public int getRowDimension() {
  119.                 return sqrtDiag.getDimension();
  120.             }

  121.             /** {@inheritDoc} */
  122.             @Override
  123.             public int getColumnDimension() {
  124.                 return sqrtDiag.getDimension();
  125.             }
  126.         };
  127.     }
  128. }