1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
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
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82 public class CMAESOptimizer
83 extends MultivariateOptimizer {
84
85
86
87
88
89
90
91
92 private int lambda;
93
94
95
96
97
98
99
100 private final boolean isActiveCMA;
101
102
103
104
105 private final int checkFeasableCount;
106
107
108
109 private double[] inputSigma;
110
111 private int dimension;
112
113
114
115
116
117
118
119
120 private int diagonalOnly;
121
122 private boolean isMinimize = true;
123
124 private final boolean generateStatistics;
125
126
127
128 private final int maxIterations;
129
130 private final double stopFitness;
131
132 private double stopTolUpX;
133
134 private double stopTolX;
135
136 private double stopTolFun;
137
138 private double stopTolHistFun;
139
140
141
142 private int mu;
143
144 private double logMu2;
145
146 private RealMatrix weights;
147
148 private double mueff;
149
150
151
152 private double sigma;
153
154 private double cc;
155
156 private double cs;
157
158 private double damps;
159
160 private double ccov1;
161
162 private double ccovmu;
163
164 private double chiN;
165
166 private double ccov1Sep;
167
168 private double ccovmuSep;
169
170
171
172 private RealMatrix xmean;
173
174 private RealMatrix pc;
175
176 private RealMatrix ps;
177
178 private double normps;
179
180 private RealMatrix B;
181
182 private RealMatrix D;
183
184 private RealMatrix BD;
185
186 private RealMatrix diagD;
187
188 private RealMatrix C;
189
190 private RealMatrix diagC;
191
192 private int iterations;
193
194
195 private double[] fitnessHistory;
196
197
198 private final RandomGenerator random;
199
200
201 private final List<Double> statisticsSigmaHistory;
202
203 private final List<RealMatrix> statisticsMeanHistory;
204
205 private final List<Double> statisticsFitnessHistory;
206
207 private final List<RealMatrix> statisticsDHistory;
208
209
210
211
212
213
214
215
216
217
218
219
220
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
246
247
248 public List<Double> getStatisticsSigmaHistory() {
249 return statisticsSigmaHistory;
250 }
251
252
253
254
255 public List<RealMatrix> getStatisticsMeanHistory() {
256 return statisticsMeanHistory;
257 }
258
259
260
261
262 public List<Double> getStatisticsFitnessHistory() {
263 return statisticsFitnessHistory;
264 }
265
266
267
268
269 public List<RealMatrix> getStatisticsDHistory() {
270 return statisticsDHistory;
271 }
272
273
274
275
276
277
278
279
280
281
282
283
284 public static class Sigma implements OptimizationData {
285
286 private final double[] s;
287
288
289
290
291
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
305
306
307 public double[] getSigma() {
308 return s.clone();
309 }
310 }
311
312
313
314
315
316
317
318
319
320
321
322 public static class PopulationSize implements OptimizationData {
323
324 private final int lambda;
325
326
327
328
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
340
341
342 public int getPopulationSize() {
343 return lambda;
344 }
345 }
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363 @Override
364 public PointValuePair optimize(OptimizationData... optData)
365 throws MathIllegalArgumentException, MathIllegalStateException {
366
367 return super.optimize(optData);
368 }
369
370
371 @Override
372 protected PointValuePair doOptimize() {
373
374 isMinimize = getGoalType().equals(GoalType.MINIMIZE);
375 final FitnessFunction fitfun = new FitnessFunction();
376 final double[] guess = getStartPoint();
377
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
389
390 generationLoop:
391 for (iterations = 1; iterations <= maxIterations; iterations++) {
392 incrementIterationCount();
393
394
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
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));
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
415 arz.setColumn(k, randn(dimension));
416 }
417 copyColumn(arxk, 0, arx, k);
418 try {
419 valuePenaltyPairs[k] = fitfun.value(arx.getColumn(k));
420 } catch (MathIllegalStateException e) {
421 break generationLoop;
422 }
423 }
424
425
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
432 final int[] arindex = sortedIndices(fitness);
433
434 final RealMatrix xold = xmean;
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
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
460
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
491 if (max(diagD) / min(diagD) > 1e7) {
492 break generationLoop;
493 }
494
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
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
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
527
528
529
530
531
532
533
534
535 @Override
536 protected void parseOptimizationData(OptimizationData... optData) {
537
538 super.parseOptimizationData(optData);
539
540
541
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
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
582
583
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
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);
597
598
599 stopTolUpX = 1e3 * max(insigma);
600 stopTolX = 1e-11 * max(insigma);
601 stopTolFun = 1e-12;
602 stopTolHistFun = 1e-13;
603
604
605 mu = lambda / 2;
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;
617
618
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;
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
634 xmean = MatrixUtils.createColumnRealMatrix(guess);
635 diagD = insigma.scalarMultiply(1 / sigma);
636 diagC = square(diagD);
637 pc = zeros(dimension, 1);
638 ps = zeros(dimension, 1);
639 normps = ps.getFrobeniusNorm();
640
641 B = eye(dimension, dimension);
642 D = ones(dimension, 1);
643 BD = times(B, repmat(diagD.transpose(), dimension, 1));
644 C = B.multiply(diag(square(D)).multiply(B.transpose()));
645 final int historySize = 10 + (int) (3 * 10 * dimension / (double) lambda);
646 fitnessHistory = new double[historySize];
647 for (int i = 0; i < historySize; i++) {
648 fitnessHistory[i] = Double.MAX_VALUE;
649 }
650 }
651
652
653
654
655
656
657
658
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
677
678
679
680
681
682 private void updateCovarianceDiagonalOnly(boolean hsig,
683 final RealMatrix bestArz) {
684
685 double oldFac = hsig ? 0 : ccov1Sep * cc * (2 - cc);
686 oldFac += 1 - ccov1Sep - ccovmuSep;
687 diagC = diagC.scalarMultiply(oldFac)
688 .add(square(pc).scalarMultiply(ccov1Sep))
689 .add(times(diagC, square(bestArz).multiply(weights))
690 .scalarMultiply(ccovmuSep));
691 diagD = sqrt(diagC);
692 if (diagonalOnly > 1 &&
693 iterations > diagonalOnly) {
694
695 diagonalOnly = 0;
696 B = eye(dimension, dimension);
697 BD = diag(diagD);
698 C = diag(diagC);
699 }
700 }
701
702
703
704
705
706
707
708
709
710
711
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);
720 final RealMatrix roneu = pc.multiplyTransposed(pc)
721 .scalarMultiply(ccov1);
722
723 double oldFac = hsig ? 0 : ccov1 * cc * (2 - cc);
724 oldFac += 1 - ccov1 - ccovmu;
725 if (isActiveCMA) {
726
727 negccov = (1 - ccovmu) * 0.25 * mueff /
728 (FastMath.pow(dimension + 2, 1.5) + 2 * mueff);
729
730
731 final double negminresidualvariance = 0.66;
732
733 final double negalphaold = 0.5;
734
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
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)
757 .add(arpos.scalarMultiply(
758 ccovmu + (1 - negalphaold) * negccov)
759 .multiply(times(repmat(weights, 1, dimension),
760 arpos.transpose())))
761 .subtract(Cneg.scalarMultiply(negccov));
762 } else {
763
764 C = C.scalarMultiply(oldFac)
765 .add(roneu)
766 .add(arpos.scalarMultiply(ccovmu)
767 .multiply(times(repmat(weights, 1, dimension),
768 arpos.transpose())));
769 }
770 }
771 updateBD(negccov);
772 }
773
774
775
776
777
778
779 private void updateBD(double negccov) {
780 if (ccov1 + ccovmu + negccov > 0 &&
781 (iterations % 1. / (ccov1 + ccovmu + negccov) / dimension / 10.) < 1) {
782
783 C = triu(C, 0).add(triu(C, 1).transpose());
784
785 final EigenDecompositionSymmetric eig = new EigenDecompositionSymmetric(C);
786 B = eig.getV();
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);
806 BD = times(B, repmat(diagD.transpose(), dimension, 1));
807 }
808 }
809
810
811
812
813
814
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
825
826
827
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
843
844
845
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
863
864
865 private static class DoubleIndex implements Comparable<DoubleIndex> {
866
867 private final double value;
868
869 private final int index;
870
871
872
873
874
875 DoubleIndex(double value, int index) {
876 this.value = value;
877 this.index = index;
878 }
879
880
881 @Override
882 public int compareTo(DoubleIndex o) {
883 return Double.compare(value, o.value);
884 }
885
886
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
902 @Override
903 public int hashCode() {
904 long bits = Double.doubleToLongBits(value);
905 return (int) ((1438542 ^ (bits >>> 32) ^ bits));
906 }
907 }
908
909
910
911 private static class ValuePenaltyPair {
912
913 private double value;
914
915 private double penalty;
916
917
918
919
920
921 ValuePenaltyPair(final double value, final double penalty) {
922 this.value = value;
923 this.penalty = penalty;
924 }
925 }
926
927
928
929
930
931
932 private class FitnessFunction {
933
934
935
936
937 private final boolean isRepairMode;
938
939
940
941 FitnessFunction() {
942 isRepairMode = true;
943 }
944
945
946
947
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
966
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
985
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
1006
1007
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
1020
1021
1022
1023
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
1037
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
1051
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
1066
1067
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
1081
1082
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
1096
1097
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
1111
1112
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
1126
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
1142
1143
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
1163
1164
1165
1166
1167
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
1178
1179
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
1191
1192
1193
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
1207
1208
1209
1210 private static RealMatrix zeros(int n, int m) {
1211 return new Array2DRowRealMatrix(n, m);
1212 }
1213
1214
1215
1216
1217
1218
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
1234
1235
1236
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
1251
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
1268
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
1285
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
1299
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
1313
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
1325
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
1337
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
1349
1350
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 }