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 package org.hipparchus.linear; 23 24 import org.hipparchus.analysis.function.Sqrt; 25 import org.hipparchus.exception.LocalizedCoreFormats; 26 import org.hipparchus.exception.MathIllegalArgumentException; 27 import org.hipparchus.util.MathArrays; 28 29 /** 30 * This class implements the standard Jacobi (diagonal) preconditioner. For a 31 * matrix A<sub>ij</sub>, this preconditioner is 32 * M = diag(1 / A<sub>11</sub>, 1 / A<sub>22</sub>, …). 33 */ 34 public class JacobiPreconditioner implements RealLinearOperator { 35 36 /** The diagonal coefficients of the preconditioner. */ 37 private final ArrayRealVector diag; 38 39 /** 40 * Creates a new instance of this class. 41 * 42 * @param diag the diagonal coefficients of the linear operator to be 43 * preconditioned 44 * @param deep {@code true} if a deep copy of the above array should be 45 * performed 46 */ 47 public JacobiPreconditioner(final double[] diag, final boolean deep) { 48 this.diag = new ArrayRealVector(diag, deep); 49 } 50 51 /** 52 * Creates a new instance of this class. This method extracts the diagonal 53 * coefficients of the specified linear operator. If {@code a} does not 54 * extend {@link AbstractRealMatrix}, then the coefficients of the 55 * underlying matrix are not accessible, coefficient extraction is made by 56 * matrix-vector products with the basis vectors (and might therefore take 57 * some time). With matrices, direct entry access is carried out. 58 * 59 * @param a the linear operator for which the preconditioner should be built 60 * @return the diagonal preconditioner made of the inverse of the diagonal 61 * coefficients of the specified linear operator 62 * @throws MathIllegalArgumentException if {@code a} is not square 63 */ 64 public static JacobiPreconditioner create(final RealLinearOperator a) 65 throws MathIllegalArgumentException { 66 final int n = a.getColumnDimension(); 67 if (a.getRowDimension() != n) { 68 throw new MathIllegalArgumentException(LocalizedCoreFormats.NON_SQUARE_OPERATOR, 69 a.getRowDimension(), n); 70 } 71 final double[] diag = new double[n]; 72 if (a instanceof AbstractRealMatrix) { 73 final AbstractRealMatrix m = (AbstractRealMatrix) a; 74 for (int i = 0; i < n; i++) { 75 diag[i] = m.getEntry(i, i); 76 } 77 } else { 78 final ArrayRealVector x = new ArrayRealVector(n); 79 for (int i = 0; i < n; i++) { 80 x.set(0.); 81 x.setEntry(i, 1.); 82 diag[i] = a.operate(x).getEntry(i); 83 } 84 } 85 return new JacobiPreconditioner(diag, false); 86 } 87 88 /** {@inheritDoc} */ 89 @Override 90 public int getColumnDimension() { 91 return diag.getDimension(); 92 } 93 94 /** {@inheritDoc} */ 95 @Override 96 public int getRowDimension() { 97 return diag.getDimension(); 98 } 99 100 /** {@inheritDoc} */ 101 @Override 102 public RealVector operate(final RealVector x) { 103 // Dimension check is carried out by ebeDivide 104 return new ArrayRealVector(MathArrays.ebeDivide(x.toArray(), 105 diag.toArray()), 106 false); 107 } 108 109 /** 110 * Returns the square root of {@code this} diagonal operator. More 111 * precisely, this method returns 112 * P = diag(1 / √A<sub>11</sub>, 1 / √A<sub>22</sub>, …). 113 * 114 * @return the square root of {@code this} preconditioner 115 */ 116 public RealLinearOperator sqrt() { 117 final RealVector sqrtDiag = diag.map(new Sqrt()); 118 return new RealLinearOperator() { 119 /** {@inheritDoc} */ 120 @Override 121 public RealVector operate(final RealVector x) { 122 return new ArrayRealVector(MathArrays.ebeDivide(x.toArray(), 123 sqrtDiag.toArray()), 124 false); 125 } 126 127 /** {@inheritDoc} */ 128 @Override 129 public int getRowDimension() { 130 return sqrtDiag.getDimension(); 131 } 132 133 /** {@inheritDoc} */ 134 @Override 135 public int getColumnDimension() { 136 return sqrtDiag.getDimension(); 137 } 138 }; 139 } 140 }