View Javadoc
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.optim.nonlinear.scalar.noderiv;
24  
25  import java.util.ArrayList;
26  import java.util.Arrays;
27  import java.util.List;
28  
29  import org.hipparchus.exception.LocalizedCoreFormats;
30  import org.hipparchus.exception.MathIllegalArgumentException;
31  import org.hipparchus.exception.MathIllegalStateException;
32  import org.hipparchus.linear.Array2DRowRealMatrix;
33  import org.hipparchus.linear.EigenDecompositionSymmetric;
34  import org.hipparchus.linear.MatrixUtils;
35  import org.hipparchus.linear.RealMatrix;
36  import org.hipparchus.optim.ConvergenceChecker;
37  import org.hipparchus.optim.OptimizationData;
38  import org.hipparchus.optim.PointValuePair;
39  import org.hipparchus.optim.nonlinear.scalar.GoalType;
40  import org.hipparchus.optim.nonlinear.scalar.MultivariateOptimizer;
41  import org.hipparchus.random.RandomGenerator;
42  import org.hipparchus.util.FastMath;
43  
44  /**
45   * An implementation of the active Covariance Matrix Adaptation Evolution Strategy (CMA-ES)
46   * for non-linear, non-convex, non-smooth, global function minimization.
47   * <p>
48   * The CMA-Evolution Strategy (CMA-ES) is a reliable stochastic optimization method
49   * which should be applied if derivative-based methods, e.g. quasi-Newton BFGS or
50   * conjugate gradient, fail due to a rugged search landscape (e.g. noise, local
51   * optima, outlier, etc.) of the objective function. Like a
52   * quasi-Newton method, the CMA-ES learns and applies a variable metric
53   * on the underlying search space. Unlike a quasi-Newton method, the
54   * CMA-ES neither estimates nor uses gradients, making it considerably more
55   * reliable in terms of finding a good, or even close to optimal, solution.
56   * <p>
57   * In general, on smooth objective functions the CMA-ES is roughly ten times
58   * slower than BFGS (counting objective function evaluations, no gradients provided).
59   * For up to \(n=10\) variables also the derivative-free simplex
60   * direct search method (Nelder and Mead) can be faster, but it is
61   * far less reliable than CMA-ES.
62   * <p>
63   * The CMA-ES is particularly well suited for non-separable
64   * and/or badly conditioned problems. To observe the advantage of CMA compared
65   * to a conventional evolution strategy, it will usually take about
66   * \(30 n\) function evaluations. On difficult problems the complete
67   * optimization (a single run) is expected to take <em>roughly</em> between
68   * \(30 n\) and \(300 n^2\)
69   * function evaluations.
70   * <p>
71   * This implementation is translated and adapted from the Matlab version
72   * of the CMA-ES algorithm as implemented in module {@code cmaes.m} version 3.51.
73   * <p>
74   * For more information, please refer to the following links:
75   * <ul>
76   *  <li><a href="http://www.lri.fr/~hansen/cmaes.m">Matlab code</a></li>
77   *  <li><a href="http://www.lri.fr/~hansen/cmaesintro.html">Introduction to CMA-ES</a></li>
78   *  <li><a href="http://en.wikipedia.org/wiki/CMA-ES">Wikipedia</a></li>
79   * </ul>
80   *
81   */
82  public class CMAESOptimizer
83      extends MultivariateOptimizer {
84      // global search parameters
85      /**
86       * Population size, offspring number. The primary strategy parameter to play
87       * with, which can be increased from its default value. Increasing the
88       * population size improves global search properties in exchange to speed.
89       * Speed decreases, as a rule, at most linearly with increasing population
90       * size. It is advisable to begin with the default small population size.
91       */
92      private int lambda; // population size
93      /**
94       * Covariance update mechanism, default is active CMA. isActiveCMA = true
95       * turns on "active CMA" with a negative update of the covariance matrix and
96       * checks for positive definiteness. OPTS.CMA.active = 2 does not check for
97       * pos. def. and is numerically faster. Active CMA usually speeds up the
98       * adaptation.
99       */
100     private final boolean isActiveCMA;
101     /**
102      * Determines how often a new random offspring is generated in case it is
103      * not feasible / beyond the defined limits, default is 0.
104      */
105     private final int checkFeasableCount;
106     /**
107      * @see Sigma
108      */
109     private double[] inputSigma;
110     /** Number of objective variables/problem dimension */
111     private int dimension;
112     /**
113      * Defines the number of initial iterations, where the covariance matrix
114      * remains diagonal and the algorithm has internally linear time complexity.
115      * diagonalOnly = 1 means keeping the covariance matrix always diagonal and
116      * this setting also exhibits linear space complexity. This can be
117      * particularly useful for dimension > 100.
118      * @see <a href="http://hal.archives-ouvertes.fr/inria-00287367/en">A Simple Modification in CMA-ES</a>
119      */
120     private int diagonalOnly;
121     /** Number of objective variables/problem dimension */
122     private boolean isMinimize = true;
123     /** Indicates whether statistic data is collected. */
124     private final boolean generateStatistics;
125 
126     // termination criteria
127     /** Maximal number of iterations allowed. */
128     private final int maxIterations;
129     /** Limit for fitness value. */
130     private final double stopFitness;
131     /** Stop if x-changes larger stopTolUpX. */
132     private double stopTolUpX;
133     /** Stop if x-change smaller stopTolX. */
134     private double stopTolX;
135     /** Stop if fun-changes smaller stopTolFun. */
136     private double stopTolFun;
137     /** Stop if back fun-changes smaller stopTolHistFun. */
138     private double stopTolHistFun;
139 
140     // selection strategy parameters
141     /** Number of parents/points for recombination. */
142     private int mu; //
143     /** log(mu + 0.5), stored for efficiency. */
144     private double logMu2;   // NOPMD - using a field here is for performance reasons
145     /** Array for weighted recombination. */
146     private RealMatrix weights;
147     /** Variance-effectiveness of sum w_i x_i. */
148     private double mueff; //
149 
150     // dynamic strategy parameters and constants
151     /** Overall standard deviation - search volume. */
152     private double sigma;
153     /** Cumulation constant. */
154     private double cc;
155     /** Cumulation constant for step-size. */
156     private double cs;
157     /** Damping for step-size. */
158     private double damps;
159     /** Learning rate for rank-one update. */
160     private double ccov1;
161     /** Learning rate for rank-mu update' */
162     private double ccovmu;
163     /** Expectation of ||N(0,I)|| == norm(randn(N,1)). */
164     private double chiN;
165     /** Learning rate for rank-one update - diagonalOnly */
166     private double ccov1Sep;
167     /** Learning rate for rank-mu update - diagonalOnly */
168     private double ccovmuSep;
169 
170     // CMA internal values - updated each generation
171     /** Objective variables. */
172     private RealMatrix xmean;
173     /** Evolution path. */
174     private RealMatrix pc;
175     /** Evolution path for sigma. */
176     private RealMatrix ps;
177     /** Norm of ps, stored for efficiency. */
178     private double normps;
179     /** Coordinate system. */
180     private RealMatrix B;
181     /** Scaling. */
182     private RealMatrix D; // NOPMD
183     /** B*D, stored for efficiency. */
184     private RealMatrix BD;
185     /** Diagonal of sqrt(D), stored for efficiency. */
186     private RealMatrix diagD;
187     /** Covariance matrix. */
188     private RealMatrix C;
189     /** Diagonal of C, used for diagonalOnly. */
190     private RealMatrix diagC;
191     /** Number of iterations already performed. */
192     private int iterations;
193 
194     /** History queue of best values. */
195     private double[] fitnessHistory;
196 
197     /** Random generator. */
198     private final RandomGenerator random;
199 
200     /** History of sigma values. */
201     private final List<Double> statisticsSigmaHistory;
202     /** History of mean matrix. */
203     private final List<RealMatrix> statisticsMeanHistory;
204     /** History of fitness values. */
205     private final List<Double> statisticsFitnessHistory;
206     /** History of D matrix. */
207     private final List<RealMatrix> statisticsDHistory;
208 
209     /** Simple constructor.
210      * @param maxIterations Maximal number of iterations.
211      * @param stopFitness Whether to stop if objective function value is smaller than
212      * {@code stopFitness}.
213      * @param isActiveCMA Chooses the covariance matrix update method.
214      * @param diagonalOnly Number of initial iterations, where the covariance matrix
215      * remains diagonal.
216      * @param checkFeasableCount Determines how often new random objective variables are
217      * generated in case they are out of bounds.
218      * @param random Random generator.
219      * @param generateStatistics Whether statistic data is collected.
220      * @param checker Convergence checker.
221      *
222      */
223     public CMAESOptimizer(int maxIterations,
224                           double stopFitness,
225                           boolean isActiveCMA,
226                           int diagonalOnly,
227                           int checkFeasableCount,
228                           RandomGenerator random,
229                           boolean generateStatistics,
230                           ConvergenceChecker<PointValuePair> checker) {
231         super(checker);
232         this.maxIterations = maxIterations;
233         this.stopFitness = stopFitness;
234         this.isActiveCMA = isActiveCMA;
235         this.diagonalOnly = diagonalOnly;
236         this.checkFeasableCount = checkFeasableCount;
237         this.random = random;
238         this.generateStatistics = generateStatistics;
239         this.statisticsSigmaHistory = new ArrayList<>();
240         this.statisticsMeanHistory = new ArrayList<>();
241         this.statisticsFitnessHistory = new ArrayList<>();
242         this.statisticsDHistory = new ArrayList<>();
243     }
244 
245     /** Get history of sigma values.
246      * @return History of sigma values
247      */
248     public List<Double> getStatisticsSigmaHistory() {
249         return statisticsSigmaHistory;
250     }
251 
252     /** Get history of mean matrix.
253      * @return History of mean matrix
254      */
255     public List<RealMatrix> getStatisticsMeanHistory() {
256         return statisticsMeanHistory;
257     }
258 
259     /** Get history of fitness values.
260      * @return History of fitness values
261      */
262     public List<Double> getStatisticsFitnessHistory() {
263         return statisticsFitnessHistory;
264     }
265 
266     /** Get history of D matrix.
267      * @return History of D matrix
268      */
269     public List<RealMatrix> getStatisticsDHistory() {
270         return statisticsDHistory;
271     }
272 
273     /**
274      * Input sigma values.
275      * They define the initial coordinate-wise standard deviations for
276      * sampling new search points around the initial guess.
277      * It is suggested to set them to the estimated distance from the
278      * initial to the desired optimum.
279      * Small values induce the search to be more local (and very small
280      * values are more likely to find a local optimum close to the initial
281      * guess).
282      * Too small values might however lead to early termination.
283      */
284     public static class Sigma implements OptimizationData {
285         /** Sigma values. */
286         private final double[] s;
287 
288         /** Simple constructor.
289          * @param s Sigma values.
290          * @throws MathIllegalArgumentException if any of the array entries is smaller
291          * than zero.
292          */
293         public Sigma(double[] s)
294             throws MathIllegalArgumentException {
295             for (double v : s) {
296                 if (v < 0) {
297                     throw new MathIllegalArgumentException(LocalizedCoreFormats.NUMBER_TOO_SMALL, v, 0);
298                 }
299             }
300 
301             this.s = s.clone();
302         }
303 
304         /** Get sigma values.
305          * @return the sigma values
306          */
307         public double[] getSigma() {
308             return s.clone();
309         }
310     }
311 
312     /**
313      * Population size.
314      * The number of offspring is the primary strategy parameter.
315      * In the absence of better clues, a good default could be an
316      * integer close to {@code 4 + 3 ln(n)}, where {@code n} is the
317      * number of optimized parameters.
318      * Increasing the population size improves global search properties
319      * at the expense of speed (which in general decreases at most
320      * linearly with increasing population size).
321      */
322     public static class PopulationSize implements OptimizationData {
323         /** Population size. */
324         private final int lambda;
325 
326         /** Simple constructor.
327          * @param size Population size.
328          * @throws MathIllegalArgumentException if {@code size <= 0}.
329          */
330         public PopulationSize(int size)
331             throws MathIllegalArgumentException {
332             if (size <= 0) {
333                 throw new MathIllegalArgumentException(LocalizedCoreFormats.NUMBER_TOO_SMALL_BOUND_EXCLUDED,
334                                                        size, 0);
335             }
336             lambda = size;
337         }
338 
339         /** Get population size.
340          * @return the population size
341          */
342         public int getPopulationSize() {
343             return lambda;
344         }
345     }
346 
347     /**
348      * {@inheritDoc}
349      *
350      * @param optData Optimization data. In addition to those documented in
351      * {@link MultivariateOptimizer#parseOptimizationData(OptimizationData[])
352      * MultivariateOptimizer}, this method will register the following data:
353      * <ul>
354      *  <li>{@link Sigma}</li>
355      *  <li>{@link PopulationSize}</li>
356      * </ul>
357      * @return {@inheritDoc}
358      * @throws MathIllegalStateException if the maximal number of
359      * evaluations is exceeded.
360      * @throws MathIllegalArgumentException if the initial guess, target, and weight
361      * arguments have inconsistent dimensions.
362      */
363     @Override
364     public PointValuePair optimize(OptimizationData... optData)
365         throws MathIllegalArgumentException, MathIllegalStateException {
366         // Set up base class and perform computation.
367         return super.optimize(optData);
368     }
369 
370     /** {@inheritDoc} */
371     @Override
372     protected PointValuePair doOptimize() {
373          // -------------------- Initialization --------------------------------
374         isMinimize = getGoalType().equals(GoalType.MINIMIZE);
375         final FitnessFunction fitfun = new FitnessFunction();
376         final double[] guess = getStartPoint();
377         // number of objective variables/problem dimension
378         dimension = guess.length;
379         initializeCMA(guess);
380         ValuePenaltyPair valuePenalty = fitfun.value(guess);
381         double bestValue = valuePenalty.value+valuePenalty.penalty;
382         push(fitnessHistory, bestValue);
383         PointValuePair optimum
384             = new PointValuePair(getStartPoint(),
385                                  isMinimize ? bestValue : -bestValue);
386         PointValuePair lastResult = null;
387 
388         // -------------------- Generation Loop --------------------------------
389 
390         generationLoop:
391         for (iterations = 1; iterations <= maxIterations; iterations++) {
392             incrementIterationCount();
393 
394             // Generate and evaluate lambda offspring
395             final RealMatrix arz = randn1(dimension, lambda);
396             final RealMatrix arx = zeros(dimension, lambda);
397             final double[] fitness = new double[lambda];
398             final ValuePenaltyPair[] valuePenaltyPairs = new ValuePenaltyPair[lambda];
399             // generate random offspring
400             for (int k = 0; k < lambda; k++) {
401                 RealMatrix arxk = null;
402                 for (int i = 0; i < checkFeasableCount + 1; i++) {
403                     if (diagonalOnly <= 0) {
404                         arxk = xmean.add(BD.multiply(arz.getColumnMatrix(k))
405                                          .scalarMultiply(sigma)); // m + sig * Normal(0,C)
406                     } else {
407                         arxk = xmean.add(times(diagD,arz.getColumnMatrix(k))
408                                          .scalarMultiply(sigma));
409                     }
410                     if (i >= checkFeasableCount ||
411                         fitfun.isFeasible(arxk.getColumn(0))) {
412                         break;
413                     }
414                     // regenerate random arguments for row
415                     arz.setColumn(k, randn(dimension));
416                 }
417                 copyColumn(arxk, 0, arx, k);
418                 try {
419                     valuePenaltyPairs[k] = fitfun.value(arx.getColumn(k)); // compute fitness
420                 } catch (MathIllegalStateException e) {
421                     break generationLoop;
422                 }
423             }
424 
425             // Compute fitnesses by adding value and penalty after scaling by value range.
426             double valueRange = valueRange(valuePenaltyPairs);
427             for (int iValue=0;iValue<valuePenaltyPairs.length;iValue++) {
428                  fitness[iValue] = valuePenaltyPairs[iValue].value + valuePenaltyPairs[iValue].penalty*valueRange;
429             }
430 
431             // Sort by fitness and compute weighted mean into xmean
432             final int[] arindex = sortedIndices(fitness);
433             // Calculate new xmean, this is selection and recombination
434             final RealMatrix xold = xmean; // for speed up of Eq. (2) and (3)
435             final RealMatrix bestArx = selectColumns(arx, Arrays.copyOf(arindex, mu));
436             xmean = bestArx.multiply(weights);
437             final RealMatrix bestArz = selectColumns(arz, Arrays.copyOf(arindex, mu));
438             final RealMatrix zmean = bestArz.multiply(weights);
439             final boolean hsig = updateEvolutionPaths(zmean, xold);
440             if (diagonalOnly <= 0) {
441                 updateCovariance(hsig, bestArx, arz, arindex, xold);
442             } else {
443                 updateCovarianceDiagonalOnly(hsig, bestArz);
444             }
445             // Adapt step size sigma - Eq. (5)
446             sigma *= FastMath.exp(FastMath.min(1, (normps/chiN - 1) * cs / damps));
447             final double bestFitness = fitness[arindex[0]];
448             final double worstFitness = fitness[arindex[arindex.length - 1]];
449             if (bestValue > bestFitness) {
450                 bestValue = bestFitness;
451                 lastResult = optimum;
452                 optimum = new PointValuePair(fitfun.repair(bestArx.getColumn(0)),
453                                              isMinimize ? bestFitness : -bestFitness);
454                 if (getConvergenceChecker() != null && lastResult != null &&
455                     getConvergenceChecker().converged(iterations, optimum, lastResult)) {
456                     break generationLoop;
457                 }
458             }
459             // handle termination criteria
460             // Break, if fitness is good enough
461             if (stopFitness != 0 && bestFitness < (isMinimize ? stopFitness : -stopFitness)) {
462                 break generationLoop;
463             }
464             final double[] sqrtDiagC = sqrt(diagC).getColumn(0);
465             final double[] pcCol = pc.getColumn(0);
466             for (int i = 0; i < dimension; i++) {
467                 if (sigma * FastMath.max(FastMath.abs(pcCol[i]), sqrtDiagC[i]) > stopTolX) {
468                     break;
469                 }
470                 if (i >= dimension - 1) {
471                     break generationLoop;
472                 }
473             }
474             for (int i = 0; i < dimension; i++) {
475                 if (sigma * sqrtDiagC[i] > stopTolUpX) {
476                     break generationLoop;
477                 }
478             }
479             final double historyBest = min(fitnessHistory);
480             final double historyWorst = max(fitnessHistory);
481             if (iterations > 2 &&
482                 FastMath.max(historyWorst, worstFitness) -
483                 FastMath.min(historyBest, bestFitness) < stopTolFun) {
484                 break generationLoop;
485             }
486             if (iterations > fitnessHistory.length &&
487                 historyWorst - historyBest < stopTolHistFun) {
488                 break generationLoop;
489             }
490             // condition number of the covariance matrix exceeds 1e14
491             if (max(diagD) / min(diagD) > 1e7) {
492                 break generationLoop;
493             }
494             // user defined termination
495             if (getConvergenceChecker() != null) {
496                 final PointValuePair current
497                     = new PointValuePair(bestArx.getColumn(0),
498                                          isMinimize ? bestFitness : -bestFitness);
499                 if (lastResult != null &&
500                     getConvergenceChecker().converged(iterations, current, lastResult)) {
501                     break generationLoop;
502                     }
503                 lastResult = current;
504             }
505             // Adjust step size in case of equal function values (flat fitness)
506             if (bestValue == fitness[arindex[(int)(0.1+lambda/4.)]]) {
507                 sigma *= FastMath.exp(0.2 + cs / damps);
508             }
509             if (iterations > 2 && FastMath.max(historyWorst, bestFitness) -
510                 FastMath.min(historyBest, bestFitness) == 0) {
511                 sigma *= FastMath.exp(0.2 + cs / damps);
512             }
513             // store best in history
514             push(fitnessHistory,bestFitness);
515             if (generateStatistics) {
516                 statisticsSigmaHistory.add(sigma);
517                 statisticsFitnessHistory.add(bestFitness);
518                 statisticsMeanHistory.add(xmean.transpose());
519                 statisticsDHistory.add(diagD.transpose().scalarMultiply(1E5));
520             }
521         }
522         return optimum;
523     }
524 
525     /**
526      * Scans the list of (required and optional) optimization data that
527      * characterize the problem.
528      *
529      * @param optData Optimization data. The following data will be looked for:
530      * <ul>
531      *  <li>{@link Sigma}</li>
532      *  <li>{@link PopulationSize}</li>
533      * </ul>
534      */
535     @Override
536     protected void parseOptimizationData(OptimizationData... optData) {
537         // Allow base class to register its own data.
538         super.parseOptimizationData(optData);
539 
540         // The existing values (as set by the previous call) are reused if
541         // not provided in the argument list.
542         for (OptimizationData data : optData) {
543             if (data instanceof Sigma) {
544                 inputSigma = ((Sigma) data).getSigma();
545                 continue;
546             }
547             if (data instanceof PopulationSize) {
548                 lambda = ((PopulationSize) data).getPopulationSize();
549                 continue;
550             }
551         }
552 
553         checkParameters();
554     }
555 
556     /**
557      * Checks dimensions and values of boundaries and inputSigma if defined.
558      */
559     private void checkParameters() {
560         if (inputSigma != null) {
561             final double[] init = getStartPoint();
562 
563             if (inputSigma.length != init.length) {
564                 throw new MathIllegalArgumentException(LocalizedCoreFormats.DIMENSIONS_MISMATCH,
565                                                        inputSigma.length, init.length);
566             }
567 
568             final double[] lB = getLowerBound();
569             final double[] uB = getUpperBound();
570 
571             for (int i = 0; i < init.length; i++) {
572                 if (inputSigma[i] > uB[i] - lB[i]) {
573                     throw new MathIllegalArgumentException(LocalizedCoreFormats.OUT_OF_RANGE_SIMPLE,
574                                                            inputSigma[i], 0, uB[i] - lB[i]);
575                 }
576             }
577         }
578     }
579 
580     /**
581      * Initialization of the dynamic search parameters
582      *
583      * @param guess Initial guess for the arguments of the fitness function.
584      */
585     private void initializeCMA(double[] guess) {
586         if (lambda <= 0) {
587             throw new MathIllegalArgumentException(LocalizedCoreFormats.NUMBER_TOO_SMALL_BOUND_EXCLUDED,
588                                                    lambda, 0);
589         }
590         // initialize sigma
591         final double[][] sigmaArray = new double[guess.length][1];
592         for (int i = 0; i < guess.length; i++) {
593             sigmaArray[i][0] = inputSigma[i];
594         }
595         final RealMatrix insigma = new Array2DRowRealMatrix(sigmaArray, false);
596         sigma = max(insigma); // overall standard deviation
597 
598         // initialize termination criteria
599         stopTolUpX = 1e3 * max(insigma);
600         stopTolX = 1e-11 * max(insigma);
601         stopTolFun = 1e-12;
602         stopTolHistFun = 1e-13;
603 
604         // initialize selection strategy parameters
605         mu = lambda / 2; // number of parents/points for recombination
606         logMu2 = FastMath.log(mu + 0.5);
607         weights = log(sequence(1, mu, 1)).scalarMultiply(-1).scalarAdd(logMu2);
608         double sumw = 0;
609         double sumwq = 0;
610         for (int i = 0; i < mu; i++) {
611             double w = weights.getEntry(i, 0);
612             sumw += w;
613             sumwq += w * w;
614         }
615         weights = weights.scalarMultiply(1 / sumw);
616         mueff = sumw * sumw / sumwq; // variance-effectiveness of sum w_i x_i
617 
618         // initialize dynamic strategy parameters and constants
619         cc = (4 + mueff / dimension) /
620                 (dimension + 4 + 2 * mueff / dimension);
621         cs = (mueff + 2) / (dimension + mueff + 3.);
622         damps = (1 + 2 * FastMath.max(0, FastMath.sqrt((mueff - 1) /
623                                                        (dimension + 1)) - 1)) *
624             FastMath.max(0.3,
625                          1 - dimension / (1e-6 + maxIterations)) + cs; // minor increment
626         ccov1 = 2 / ((dimension + 1.3) * (dimension + 1.3) + mueff);
627         ccovmu = FastMath.min(1 - ccov1, 2 * (mueff - 2 + 1 / mueff) /
628                               ((dimension + 2) * (dimension + 2) + mueff));
629         ccov1Sep = FastMath.min(1, ccov1 * (dimension + 1.5) / 3);
630         ccovmuSep = FastMath.min(1 - ccov1, ccovmu * (dimension + 1.5) / 3);
631         chiN = FastMath.sqrt(dimension) *
632                 (1 - 1 / ((double) 4 * dimension) + 1 / ((double) 21 * dimension * dimension));
633         // intialize CMA internal values - updated each generation
634         xmean = MatrixUtils.createColumnRealMatrix(guess); // objective variables
635         diagD = insigma.scalarMultiply(1 / sigma);
636         diagC = square(diagD);
637         pc = zeros(dimension, 1); // evolution paths for C and sigma
638         ps = zeros(dimension, 1); // B defines the coordinate system
639         normps = ps.getFrobeniusNorm();
640 
641         B = eye(dimension, dimension);
642         D = ones(dimension, 1); // diagonal D defines the scaling
643         BD = times(B, repmat(diagD.transpose(), dimension, 1));
644         C = B.multiply(diag(square(D)).multiply(B.transpose())); // covariance
645         final int historySize = 10 + (int) (3 * 10 * dimension / (double) lambda);
646         fitnessHistory = new double[historySize]; // history of fitness values
647         for (int i = 0; i < historySize; i++) {
648             fitnessHistory[i] = Double.MAX_VALUE;
649         }
650     }
651 
652     /**
653      * Update of the evolution paths ps and pc.
654      *
655      * @param zmean Weighted row matrix of the gaussian random numbers generating
656      * the current offspring.
657      * @param xold xmean matrix of the previous generation.
658      * @return hsig flag indicating a small correction.
659      */
660     private boolean updateEvolutionPaths(RealMatrix zmean, RealMatrix xold) {
661         ps = ps.scalarMultiply(1 - cs).add(
662                 B.multiply(zmean).scalarMultiply(
663                         FastMath.sqrt(cs * (2 - cs) * mueff)));
664         normps = ps.getFrobeniusNorm();
665         final boolean hsig = normps /
666             FastMath.sqrt(1 - FastMath.pow(1 - cs, 2 * iterations)) /
667             chiN < 1.4 + 2 / ((double) dimension + 1);
668         pc = pc.scalarMultiply(1 - cc);
669         if (hsig) {
670             pc = pc.add(xmean.subtract(xold).scalarMultiply(FastMath.sqrt(cc * (2 - cc) * mueff) / sigma));
671         }
672         return hsig;
673     }
674 
675     /**
676      * Update of the covariance matrix C for diagonalOnly > 0
677      *
678      * @param hsig Flag indicating a small correction.
679      * @param bestArz Fitness-sorted matrix of the gaussian random values of the
680      * current offspring.
681      */
682     private void updateCovarianceDiagonalOnly(boolean hsig,
683                                               final RealMatrix bestArz) {
684         // minor correction if hsig==false
685         double oldFac = hsig ? 0 : ccov1Sep * cc * (2 - cc);
686         oldFac += 1 - ccov1Sep - ccovmuSep;
687         diagC = diagC.scalarMultiply(oldFac) // regard old matrix
688             .add(square(pc).scalarMultiply(ccov1Sep)) // plus rank one update
689             .add(times(diagC, square(bestArz).multiply(weights)) // plus rank mu update
690                  .scalarMultiply(ccovmuSep));
691         diagD = sqrt(diagC); // replaces eig(C)
692         if (diagonalOnly > 1 &&
693             iterations > diagonalOnly) {
694             // full covariance matrix from now on
695             diagonalOnly = 0;
696             B = eye(dimension, dimension);
697             BD = diag(diagD);
698             C = diag(diagC);
699         }
700     }
701 
702     /**
703      * Update of the covariance matrix C.
704      *
705      * @param hsig Flag indicating a small correction.
706      * @param bestArx Fitness-sorted matrix of the argument vectors producing the
707      * current offspring.
708      * @param arz Unsorted matrix containing the gaussian random values of the
709      * current offspring.
710      * @param arindex Indices indicating the fitness-order of the current offspring.
711      * @param xold xmean matrix of the previous generation.
712      */
713     private void updateCovariance(boolean hsig, final RealMatrix bestArx,
714                                   final RealMatrix arz, final int[] arindex,
715                                   final RealMatrix xold) {
716         double negccov = 0;
717         if (ccov1 + ccovmu > 0) {
718             final RealMatrix arpos = bestArx.subtract(repmat(xold, 1, mu))
719                 .scalarMultiply(1 / sigma); // mu difference vectors
720             final RealMatrix roneu = pc.multiplyTransposed(pc)
721                 .scalarMultiply(ccov1); // rank one update
722             // minor correction if hsig==false
723             double oldFac = hsig ? 0 : ccov1 * cc * (2 - cc);
724             oldFac += 1 - ccov1 - ccovmu;
725             if (isActiveCMA) {
726                 // Adapt covariance matrix C active CMA
727                 negccov = (1 - ccovmu) * 0.25 * mueff /
728                     (FastMath.pow(dimension + 2, 1.5) + 2 * mueff);
729                 // keep at least 0.66 in all directions, small popsize are most
730                 // critical
731                 final double negminresidualvariance = 0.66;
732                 // where to make up for the variance loss
733                 final double negalphaold = 0.5;
734                 // prepare vectors, compute negative updating matrix Cneg
735                 final int[] arReverseIndex = reverse(arindex);
736                 RealMatrix arzneg = selectColumns(arz, Arrays.copyOf(arReverseIndex, mu));
737                 RealMatrix arnorms = sqrt(sumRows(square(arzneg)));
738                 final int[] idxnorms = sortedIndices(arnorms.getRow(0));
739                 final RealMatrix arnormsSorted = selectColumns(arnorms, idxnorms);
740                 final int[] idxReverse = reverse(idxnorms);
741                 final RealMatrix arnormsReverse = selectColumns(arnorms, idxReverse);
742                 arnorms = divide(arnormsReverse, arnormsSorted);
743                 final int[] idxInv = inverse(idxnorms);
744                 final RealMatrix arnormsInv = selectColumns(arnorms, idxInv);
745                 // check and set learning rate negccov
746                 final double negcovMax = (1 - negminresidualvariance) /
747                     square(arnormsInv).multiply(weights).getEntry(0, 0);
748                 if (negccov > negcovMax) {
749                     negccov = negcovMax;
750                 }
751                 arzneg = times(arzneg, repmat(arnormsInv, dimension, 1));
752                 final RealMatrix artmp = BD.multiply(arzneg);
753                 final RealMatrix Cneg = artmp.multiply(diag(weights)).multiply(artmp.transpose());
754                 oldFac += negalphaold * negccov;
755                 C = C.scalarMultiply(oldFac)
756                     .add(roneu) // regard old matrix
757                     .add(arpos.scalarMultiply( // plus rank one update
758                                               ccovmu + (1 - negalphaold) * negccov) // plus rank mu update
759                          .multiply(times(repmat(weights, 1, dimension),
760                                          arpos.transpose())))
761                     .subtract(Cneg.scalarMultiply(negccov));
762             } else {
763                 // Adapt covariance matrix C - nonactive
764                 C = C.scalarMultiply(oldFac) // regard old matrix
765                     .add(roneu) // plus rank one update
766                     .add(arpos.scalarMultiply(ccovmu) // plus rank mu update
767                          .multiply(times(repmat(weights, 1, dimension),
768                                          arpos.transpose())));
769             }
770         }
771         updateBD(negccov);
772     }
773 
774     /**
775      * Update B and D from C.
776      *
777      * @param negccov Negative covariance factor.
778      */
779     private void updateBD(double negccov) {
780         if (ccov1 + ccovmu + negccov > 0 &&
781             (iterations % 1. / (ccov1 + ccovmu + negccov) / dimension / 10.) < 1) {
782             // to achieve O(N^2)
783             C = triu(C, 0).add(triu(C, 1).transpose());
784             // enforce symmetry to prevent complex numbers
785             final EigenDecompositionSymmetric eig = new EigenDecompositionSymmetric(C);
786             B = eig.getV(); // eigen decomposition, B==normalized eigenvectors
787             D = eig.getD();
788             diagD = diag(D);
789             if (min(diagD) <= 0) {
790                 for (int i = 0; i < dimension; i++) {
791                     if (diagD.getEntry(i, 0) < 0) {
792                         diagD.setEntry(i, 0, 0);
793                     }
794                 }
795                 final double tfac = max(diagD) / 1e14;
796                 C = C.add(eye(dimension, dimension).scalarMultiply(tfac));
797                 diagD = diagD.add(ones(dimension, 1).scalarMultiply(tfac));
798             }
799             if (max(diagD) > 1e14 * min(diagD)) {
800                 final double tfac = max(diagD) / 1e14 - min(diagD);
801                 C = C.add(eye(dimension, dimension).scalarMultiply(tfac));
802                 diagD = diagD.add(ones(dimension, 1).scalarMultiply(tfac));
803             }
804             diagC = diag(C);
805             diagD = sqrt(diagD); // D contains standard deviations now
806             BD = times(B, repmat(diagD.transpose(), dimension, 1)); // O(n^2)
807         }
808     }
809 
810     /**
811      * Pushes the current best fitness value in a history queue.
812      *
813      * @param vals History queue.
814      * @param val Current best fitness value.
815      */
816     private static void push(double[] vals, double val) {
817         for (int i = vals.length-1; i > 0; i--) {
818             vals[i] = vals[i-1];
819         }
820         vals[0] = val;
821     }
822 
823     /**
824      * Sorts fitness values.
825      *
826      * @param doubles Array of values to be sorted.
827      * @return a sorted array of indices pointing into doubles.
828      */
829     private int[] sortedIndices(final double[] doubles) {
830         final DoubleIndex[] dis = new DoubleIndex[doubles.length];
831         for (int i = 0; i < doubles.length; i++) {
832             dis[i] = new DoubleIndex(doubles[i], i);
833         }
834         Arrays.sort(dis);
835         final int[] indices = new int[doubles.length];
836         for (int i = 0; i < doubles.length; i++) {
837             indices[i] = dis[i].index;
838         }
839         return indices;
840     }
841    /**
842      * Get range of values.
843      *
844      * @param vpPairs Array of valuePenaltyPairs to get range from.
845      * @return a double equal to maximum value minus minimum value.
846      */
847     private double valueRange(final ValuePenaltyPair[] vpPairs) {
848         double max = Double.NEGATIVE_INFINITY;
849         double min = Double.MAX_VALUE;
850         for (ValuePenaltyPair vpPair:vpPairs) {
851             if (vpPair.value > max) {
852                 max = vpPair.value;
853             }
854             if (vpPair.value < min) {
855                 min = vpPair.value;
856             }
857         }
858         return max-min;
859     }
860 
861     /**
862      * Used to sort fitness values. Sorting is always in lower value first
863      * order.
864      */
865     private static class DoubleIndex implements Comparable<DoubleIndex> {
866         /** Value to compare. */
867         private final double value;
868         /** Index into sorted array. */
869         private final int index;
870 
871         /**
872          * @param value Value to compare.
873          * @param index Index into sorted array.
874          */
875         DoubleIndex(double value, int index) {
876             this.value = value;
877             this.index = index;
878         }
879 
880         /** {@inheritDoc} */
881         @Override
882         public int compareTo(DoubleIndex o) {
883             return Double.compare(value, o.value);
884         }
885 
886         /** {@inheritDoc} */
887         @Override
888         public boolean equals(Object other) {
889 
890             if (this == other) {
891                 return true;
892             }
893 
894             if (other instanceof DoubleIndex) {
895                 return Double.compare(value, ((DoubleIndex) other).value) == 0;
896             }
897 
898             return false;
899         }
900 
901         /** {@inheritDoc} */
902         @Override
903         public int hashCode() {
904             long bits = Double.doubleToLongBits(value);
905             return (int) ((1438542 ^ (bits >>> 32) ^ bits));
906         }
907     }
908     /**
909      * Stores the value and penalty (for repair of out of bounds point).
910      */
911     private static class ValuePenaltyPair {
912         /** Objective function value. */
913         private double value;
914         /** Penalty value for repair of out out of bounds points. */
915         private double penalty;
916 
917         /**
918          * @param value Function value.
919          * @param penalty Out-of-bounds penalty.
920         */
921         ValuePenaltyPair(final double value, final double penalty) {
922             this.value   = value;
923             this.penalty = penalty;
924         }
925     }
926 
927 
928     /**
929      * Normalizes fitness values to the range [0,1]. Adds a penalty to the
930      * fitness value if out of range.
931      */
932     private class FitnessFunction {
933         /**
934          * Flag indicating whether the objective variables are forced into their
935          * bounds if defined
936          */
937         private final boolean isRepairMode;
938 
939         /** Simple constructor.
940          */
941         FitnessFunction() {
942             isRepairMode = true;
943         }
944 
945         /**
946          * @param point Normalized objective variables.
947          * @return the objective value + penalty for violated bounds.
948          */
949         public ValuePenaltyPair value(final double[] point) {
950             double value;
951             double penalty=0.0;
952             if (isRepairMode) {
953                 double[] repaired = repair(point);
954                 value = CMAESOptimizer.this.computeObjectiveValue(repaired);
955                 penalty =  penalty(point, repaired);
956             } else {
957                 value = CMAESOptimizer.this.computeObjectiveValue(point);
958             }
959             value = isMinimize ? value : -value;
960             penalty = isMinimize ? penalty : -penalty;
961             return new ValuePenaltyPair(value,penalty);
962         }
963 
964         /**
965          * @param x Normalized objective variables.
966          * @return {@code true} if in bounds.
967          */
968         public boolean isFeasible(final double[] x) {
969             final double[] lB = CMAESOptimizer.this.getLowerBound();
970             final double[] uB = CMAESOptimizer.this.getUpperBound();
971 
972             for (int i = 0; i < x.length; i++) {
973                 if (x[i] < lB[i]) {
974                     return false;
975                 }
976                 if (x[i] > uB[i]) {
977                     return false;
978                 }
979             }
980             return true;
981         }
982 
983         /**
984          * @param x Normalized objective variables.
985          * @return the repaired (i.e. all in bounds) objective variables.
986          */
987         private double[] repair(final double[] x) {
988             final double[] lB = CMAESOptimizer.this.getLowerBound();
989             final double[] uB = CMAESOptimizer.this.getUpperBound();
990 
991             final double[] repaired = new double[x.length];
992             for (int i = 0; i < x.length; i++) {
993                 if (x[i] < lB[i]) {
994                     repaired[i] = lB[i];
995                 } else if (x[i] > uB[i]) {
996                     repaired[i] = uB[i];
997                 } else {
998                     repaired[i] = x[i];
999                 }
1000             }
1001             return repaired;
1002         }
1003 
1004         /**
1005          * @param x Normalized objective variables.
1006          * @param repaired Repaired objective variables.
1007          * @return Penalty value according to the violation of the bounds.
1008          */
1009         private double penalty(final double[] x, final double[] repaired) {
1010             double penalty = 0;
1011             for (int i = 0; i < x.length; i++) {
1012                 double diff = FastMath.abs(x[i] - repaired[i]);
1013                 penalty += diff;
1014             }
1015             return isMinimize ? penalty : -penalty;
1016         }
1017     }
1018 
1019     // -----Matrix utility functions similar to the Matlab build in functions------
1020 
1021     /**
1022      * @param m Input matrix
1023      * @return Matrix representing the element-wise logarithm of m.
1024      */
1025     private static RealMatrix log(final RealMatrix m) {
1026         final double[][] d = new double[m.getRowDimension()][m.getColumnDimension()];
1027         for (int r = 0; r < m.getRowDimension(); r++) {
1028             for (int c = 0; c < m.getColumnDimension(); c++) {
1029                 d[r][c] = FastMath.log(m.getEntry(r, c));
1030             }
1031         }
1032         return new Array2DRowRealMatrix(d, false);
1033     }
1034 
1035     /**
1036      * @param m Input matrix.
1037      * @return Matrix representing the element-wise square root of m.
1038      */
1039     private static RealMatrix sqrt(final RealMatrix m) {
1040         final double[][] d = new double[m.getRowDimension()][m.getColumnDimension()];
1041         for (int r = 0; r < m.getRowDimension(); r++) {
1042             for (int c = 0; c < m.getColumnDimension(); c++) {
1043                 d[r][c] = FastMath.sqrt(m.getEntry(r, c));
1044             }
1045         }
1046         return new Array2DRowRealMatrix(d, false);
1047     }
1048 
1049     /**
1050      * @param m Input matrix.
1051      * @return Matrix representing the element-wise square of m.
1052      */
1053     private static RealMatrix square(final RealMatrix m) {
1054         final double[][] d = new double[m.getRowDimension()][m.getColumnDimension()];
1055         for (int r = 0; r < m.getRowDimension(); r++) {
1056             for (int c = 0; c < m.getColumnDimension(); c++) {
1057                 double e = m.getEntry(r, c);
1058                 d[r][c] = e * e;
1059             }
1060         }
1061         return new Array2DRowRealMatrix(d, false);
1062     }
1063 
1064     /**
1065      * @param m Input matrix 1.
1066      * @param n Input matrix 2.
1067      * @return the matrix where the elements of m and n are element-wise multiplied.
1068      */
1069     private static RealMatrix times(final RealMatrix m, final RealMatrix n) {
1070         final double[][] d = new double[m.getRowDimension()][m.getColumnDimension()];
1071         for (int r = 0; r < m.getRowDimension(); r++) {
1072             for (int c = 0; c < m.getColumnDimension(); c++) {
1073                 d[r][c] = m.getEntry(r, c) * n.getEntry(r, c);
1074             }
1075         }
1076         return new Array2DRowRealMatrix(d, false);
1077     }
1078 
1079     /**
1080      * @param m Input matrix 1.
1081      * @param n Input matrix 2.
1082      * @return Matrix where the elements of m and n are element-wise divided.
1083      */
1084     private static RealMatrix divide(final RealMatrix m, final RealMatrix n) {
1085         final double[][] d = new double[m.getRowDimension()][m.getColumnDimension()];
1086         for (int r = 0; r < m.getRowDimension(); r++) {
1087             for (int c = 0; c < m.getColumnDimension(); c++) {
1088                 d[r][c] = m.getEntry(r, c) / n.getEntry(r, c);
1089             }
1090         }
1091         return new Array2DRowRealMatrix(d, false);
1092     }
1093 
1094     /**
1095      * @param m Input matrix.
1096      * @param cols Columns to select.
1097      * @return Matrix representing the selected columns.
1098      */
1099     private static RealMatrix selectColumns(final RealMatrix m, final int[] cols) {
1100         final double[][] d = new double[m.getRowDimension()][cols.length];
1101         for (int r = 0; r < m.getRowDimension(); r++) {
1102             for (int c = 0; c < cols.length; c++) {
1103                 d[r][c] = m.getEntry(r, cols[c]);
1104             }
1105         }
1106         return new Array2DRowRealMatrix(d, false);
1107     }
1108 
1109     /**
1110      * @param m Input matrix.
1111      * @param k Diagonal position.
1112      * @return Upper triangular part of matrix.
1113      */
1114     private static RealMatrix triu(final RealMatrix m, int k) {
1115         final double[][] d = new double[m.getRowDimension()][m.getColumnDimension()];
1116         for (int r = 0; r < m.getRowDimension(); r++) {
1117             for (int c = 0; c < m.getColumnDimension(); c++) {
1118                 d[r][c] = r <= c - k ? m.getEntry(r, c) : 0;
1119             }
1120         }
1121         return new Array2DRowRealMatrix(d, false);
1122     }
1123 
1124     /**
1125      * @param m Input matrix.
1126      * @return Row matrix representing the sums of the rows.
1127      */
1128     private static RealMatrix sumRows(final RealMatrix m) {
1129         final double[][] d = new double[1][m.getColumnDimension()];
1130         for (int c = 0; c < m.getColumnDimension(); c++) {
1131             double sum = 0;
1132             for (int r = 0; r < m.getRowDimension(); r++) {
1133                 sum += m.getEntry(r, c);
1134             }
1135             d[0][c] = sum;
1136         }
1137         return new Array2DRowRealMatrix(d, false);
1138     }
1139 
1140     /**
1141      * @param m Input matrix.
1142      * @return the diagonal n-by-n matrix if m is a column matrix or the column
1143      * matrix representing the diagonal if m is a n-by-n matrix.
1144      */
1145     private static RealMatrix diag(final RealMatrix m) {
1146         if (m.getColumnDimension() == 1) {
1147             final double[][] d = new double[m.getRowDimension()][m.getRowDimension()];
1148             for (int i = 0; i < m.getRowDimension(); i++) {
1149                 d[i][i] = m.getEntry(i, 0);
1150             }
1151             return new Array2DRowRealMatrix(d, false);
1152         } else {
1153             final double[][] d = new double[m.getRowDimension()][1];
1154             for (int i = 0; i < m.getColumnDimension(); i++) {
1155                 d[i][0] = m.getEntry(i, i);
1156             }
1157             return new Array2DRowRealMatrix(d, false);
1158         }
1159     }
1160 
1161     /**
1162      * Copies a column from m1 to m2.
1163      *
1164      * @param m1 Source matrix.
1165      * @param col1 Source column.
1166      * @param m2 Target matrix.
1167      * @param col2 Target column.
1168      */
1169     private static void copyColumn(final RealMatrix m1, int col1,
1170                                    RealMatrix m2, int col2) {
1171         for (int i = 0; i < m1.getRowDimension(); i++) {
1172             m2.setEntry(i, col2, m1.getEntry(i, col1));
1173         }
1174     }
1175 
1176     /**
1177      * @param n Number of rows.
1178      * @param m Number of columns.
1179      * @return n-by-m matrix filled with 1.
1180      */
1181     private static RealMatrix ones(int n, int m) {
1182         final double[][] d = new double[n][m];
1183         for (int r = 0; r < n; r++) {
1184             Arrays.fill(d[r], 1);
1185         }
1186         return new Array2DRowRealMatrix(d, false);
1187     }
1188 
1189     /**
1190      * @param n Number of rows.
1191      * @param m Number of columns.
1192      * @return n-by-m matrix of 0 values out of diagonal, and 1 values on
1193      * the diagonal.
1194      */
1195     private static RealMatrix eye(int n, int m) {
1196         final double[][] d = new double[n][m];
1197         for (int r = 0; r < n; r++) {
1198             if (r < m) {
1199                 d[r][r] = 1;
1200             }
1201         }
1202         return new Array2DRowRealMatrix(d, false);
1203     }
1204 
1205     /**
1206      * @param n Number of rows.
1207      * @param m Number of columns.
1208      * @return n-by-m matrix of zero values.
1209      */
1210     private static RealMatrix zeros(int n, int m) {
1211         return new Array2DRowRealMatrix(n, m);
1212     }
1213 
1214     /**
1215      * @param mat Input matrix.
1216      * @param n Number of row replicates.
1217      * @param m Number of column replicates.
1218      * @return a matrix which replicates the input matrix in both directions.
1219      */
1220     private static RealMatrix repmat(final RealMatrix mat, int n, int m) {
1221         final int rd = mat.getRowDimension();
1222         final int cd = mat.getColumnDimension();
1223         final double[][] d = new double[n * rd][m * cd];
1224         for (int r = 0; r < n * rd; r++) {
1225             for (int c = 0; c < m * cd; c++) {
1226                 d[r][c] = mat.getEntry(r % rd, c % cd);
1227             }
1228         }
1229         return new Array2DRowRealMatrix(d, false);
1230     }
1231 
1232     /**
1233      * @param start Start value.
1234      * @param end End value.
1235      * @param step Step size.
1236      * @return a sequence as column matrix.
1237      */
1238     private static RealMatrix sequence(double start, double end, double step) {
1239         final int size = (int) ((end - start) / step + 1);
1240         final double[][] d = new double[size][1];
1241         double value = start;
1242         for (int r = 0; r < size; r++) {
1243             d[r][0] = value;
1244             value += step;
1245         }
1246         return new Array2DRowRealMatrix(d, false);
1247     }
1248 
1249     /**
1250      * @param m Input matrix.
1251      * @return the maximum of the matrix element values.
1252      */
1253     private static double max(final RealMatrix m) {
1254         double max = -Double.MAX_VALUE;
1255         for (int r = 0; r < m.getRowDimension(); r++) {
1256             for (int c = 0; c < m.getColumnDimension(); c++) {
1257                 double e = m.getEntry(r, c);
1258                 if (max < e) {
1259                     max = e;
1260                 }
1261             }
1262         }
1263         return max;
1264     }
1265 
1266     /**
1267      * @param m Input matrix.
1268      * @return the minimum of the matrix element values.
1269      */
1270     private static double min(final RealMatrix m) {
1271         double min = Double.MAX_VALUE;
1272         for (int r = 0; r < m.getRowDimension(); r++) {
1273             for (int c = 0; c < m.getColumnDimension(); c++) {
1274                 double e = m.getEntry(r, c);
1275                 if (min > e) {
1276                     min = e;
1277                 }
1278             }
1279         }
1280         return min;
1281     }
1282 
1283     /**
1284      * @param m Input array.
1285      * @return the maximum of the array values.
1286      */
1287     private static double max(final double[] m) {
1288         double max = -Double.MAX_VALUE;
1289         for (double v : m) {
1290             if (max < v) {
1291                 max = v;
1292             }
1293         }
1294         return max;
1295     }
1296 
1297     /**
1298      * @param m Input array.
1299      * @return the minimum of the array values.
1300      */
1301     private static double min(final double[] m) {
1302         double min = Double.MAX_VALUE;
1303         for (double v : m) {
1304             if (min > v) {
1305                 min = v;
1306             }
1307         }
1308         return min;
1309     }
1310 
1311     /**
1312      * @param indices Input index array.
1313      * @return the inverse of the mapping defined by indices.
1314      */
1315     private static int[] inverse(final int[] indices) {
1316         final int[] inverse = new int[indices.length];
1317         for (int i = 0; i < indices.length; i++) {
1318             inverse[indices[i]] = i;
1319         }
1320         return inverse;
1321     }
1322 
1323     /**
1324      * @param indices Input index array.
1325      * @return the indices in inverse order (last is first).
1326      */
1327     private static int[] reverse(final int[] indices) {
1328         final int[] reverse = new int[indices.length];
1329         for (int i = 0; i < indices.length; i++) {
1330             reverse[i] = indices[indices.length - i - 1];
1331         }
1332         return reverse;
1333     }
1334 
1335     /**
1336      * @param size Length of random array.
1337      * @return an array of Gaussian random numbers.
1338      */
1339     private double[] randn(int size) {
1340         final double[] randn = new double[size];
1341         for (int i = 0; i < size; i++) {
1342             randn[i] = random.nextGaussian();
1343         }
1344         return randn;
1345     }
1346 
1347     /**
1348      * @param size Number of rows.
1349      * @param popSize Population size.
1350      * @return a 2-dimensional matrix of Gaussian random numbers.
1351      */
1352     private RealMatrix randn1(int size, int popSize) {
1353         final double[][] d = new double[size][popSize];
1354         for (int r = 0; r < size; r++) {
1355             for (int c = 0; c < popSize; c++) {
1356                 d[r][c] = random.nextGaussian();
1357             }
1358         }
1359         return new Array2DRowRealMatrix(d, false);
1360     }
1361 }