gaussAdapt.m 28.3 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 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 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046
function [xmin,fmin,counteval,out] = gaussAdapt(fitfun,xstart,inopts)
%------------------------------------------------------------------------
% Implementation of the Gaussian Adaptation algorithm
% When using this code please cite:
%
%          
% Christian L. Mueller
% MOSAIC group, ETH Zurich, Switzerland
%-------------------------------------------------------------------------

% dimension of the problem
N = length(xstart);

% Options defaults: Stopping criteria % (value of stop flag)
defopts.StopFitness   = -Inf;    % stop if f(xmin) < stopfitness, minimization';
defopts.MaxIter       = 1e4*(N); % maximal number of iterations/function evaluations';
defopts.TolX          = 1e-12;   % restart if history of xvals smaller TolX';
defopts.TolFun        = 1e-9;    % restart if history of funvals smaller TolFun';
defopts.TolR          = 1e-9;    % restart if step size smaller TolR';
defopts.TolCon        = 1e-9;    % restart if threshold and fitness converge';
defopts.BoundActive   = 0;       % Flag for existence of bounds
defopts.BoundPenalty  = 0;       % Flag for use of penalty term outside bounds;
defopts.LBounds       = -Inf;    % lower bounds, scalar or Nx1-vector';
defopts.UBounds       = Inf;     % upper bounds, scalar or Nx1-vector';
defopts.bRestart      = 1;       % Flag for restart activation;
defopts.ThreshRank    = 0;       % Flag for threshold based on ranks (Experimental)
defopts.PopMode       = 0;       % Flag for population mode (ToDo);
defopts.Display       = 'off';   % Display 2D landscape while running';
defopts.Plotting      = 'on';    % Plot progress while running';
defopts.VerboseModulo = 1e3;     % >=0, command line messages after every i-th iteration';
defopts.SavingModulo  = 1e2;     % >=0, saving after every i-th iteration';
defopts.bSaving       = 'on';    % [on|off] save data to file';
defopts.bSaveCov      = 0;       % save covariance matrices' ('1' results in huge files);
defopts.funArgs       = [];      % give function arguments (scalar or vectors) if necessary

% Default options for algorithmic parameters
defopts.valP    = 1/exp(1);         % 0.37... Hitting probability
defopts.r       = 1;                % Step size of the initial covariance
defopts.MaxR    = 10;               % maximal allowed step size
defopts.N_mu    = exp(1)*N;         % Mean adaptation weight
defopts.N_C     = (N+1)^2/log(N+1); % Matrix adaptation weight
defopts.N_T     = exp(1)*N;         % Constraint adaptation weight
defopts.inc_T   = 2;                % Factor for N_T increase at restart
defopts.beta    = 1/defopts.N_C;    % Step size increase/decrease factor
defopts.ss      = 1 + defopts.beta*(1-defopts.valP); % Expansion upon success
defopts.sf      = 1 - defopts.beta*(defopts.valP);   % Contraction otherwise

% Options for optimization/design-centering setting
defopts.Optim = 1;
defopts.c_T   = 1;

% ---------------------- Handling Input Parameters ----------------------

if isempty(fitfun)
    error('Objective function not determined');
end
if ~ischar(fitfun)
    error('first argument FUN must be a string');
end

if nargin < 2
    xstart = [];
end

if isempty(xstart)
    error('Initial search point, and problem dimension, not determined');
end

% Compose options opts
if nargin < 3 || isempty(inopts) % no input options available
    inopts = [];
    opts = defopts;
else
    opts = getoptions(inopts, defopts);
end


% ---------------------- Setup algorithmic Parameters ----------------------

% Options for algorithmic parameters
P = opts.valP;

if ~isfield(inopts,'N_T') && isfield(inopts,'N_C')
    opts.N_T = opts.N_C/2;
end

if ~isfield(inopts,'beta') && isfield(inopts,'N_C')
    opts.beta = 1/(opts.N_C/2);
end

if ~isfield(inopts,'ss')
    opts.ss = 1 + opts.beta*(1-opts.valP); % Expansion
end

if ~isfield(inopts,'sf')
    opts.sf = 1 - opts.beta*(opts.valP); % Contraction
end

if isfield(inopts,'UBounds') && isfield(inopts,'LBounds')
    opts.BoundActive  = 1;
    if ~isfield(inopts,'MaxR')
        opts.MaxR = 0.5*max(opts.UBounds-opts.LBounds);
    end
end


N_mu = opts.N_mu;
N_C = opts.N_C;
N_T = opts.N_T;
beta = opts.beta;
ss = opts.ss;
sf = opts.sf;
r = opts.r;
r_init = r;
inc_T = opts.inc_T;
rMax = opts.MaxR;

P_emp = 0;

% Display the current options
% opts
% pause

% Initialize dynamic (internal) strategy parameters
Q = eye(N,N);
C =  r^2*eye(N,N);

% compute function for plotting
if strcmp(opts.Plotting,'on')
    
    figure(42)
    hold on
    grid on
    lastInd = 1;
    lastSave = 1;
    lastEigVals = ones(N,1);
    
end

if strcmp(opts.Display,'on')
    
    if strcmp(fitfun,'benchmark_func')
        global initial_flag;
        initial_flag=0;
    end
    
    figure(23)
    
    x = linspace(opts.LBounds(1),opts.UBounds(1),100);
    y = linspace(opts.LBounds(2),opts.UBounds(2),100);
    f = zeros(length(x),length(y));
    
    i=1;
    j=1;
    for j=1:length(x)
        for i=1:length(y)
            if isempty(opts.funArgs)
                f(i,j)=feval(fitfun, [x(j);y(i)]);
            else
                f(i,j)=feval(fitfun, [x(j);y(i)],opts.funArgs);
            end
        end
    end
    
    [X,Y]=meshgrid(x,y);
    % meshc(Y,X,f)
    % hold on
    if strcmp(fitfun,'fnoisysphere')
        contour(X,Y,f,10)
    else
        contour(X,Y,f,100)
    end
    
    if strcmp(fitfun,'benchmark_func')
        initial_flag=0;
    end
    
    hold on
    
end

% ----------------- Setup initial settings ---------------------------

% Determine start fitness value and threshold
if isempty(opts.funArgs)
    c_T=feval(fitfun, xstart);
else
    c_T=feval(fitfun, xstart,opts.funArgs);
end

arfitness=c_T;

if (opts.Optim==0)
    if isfield(inopts,'c_T')
        c_T = opts.c_T;
        threshold = c_T;
    end
end

% For testing issues save the initial values
arfitness_init = arfitness;
c_T_init = c_T;
N_mu_init = N_mu;
N_C_init = N_C;
N_T_init = N_T;

% Number of iterations
counteval = 1;

% Number of accepted points
numAcc = 0;

% Number of accepted points for a restart
rNumAcc = 0;

% ----------------- Setup output Parameters ---------------------------

mu = xstart;
mu_old=mu;
arx = xstart;

% Restart specific bestever
currxBest = xstart;
currfBest = arfitness;

% Overall bestever
xmin = xstart;
fmin = arfitness;

% History for rank-based selection
histLen=10*N;%2*ceil(4+3*log(N));
histFAcc = zeros(histLen,1);
histXAcc = zeros(histLen,N);
histFAcc(:) = currfBest;
sortedHist = histFAcc;

if strcmp(opts.bSaving,'on')
    
    % Trace of raw function values
    xRaw=zeros(opts.MaxIter/opts.SavingModulo,N);
    xRaw(1,:)=xstart;
    fRaw=zeros(opts.MaxIter/opts.SavingModulo,1);
    fRaw(1) = c_T;
    
    xAcc=zeros(opts.MaxIter,N);
    xAcc(1,:)=xstart';
    fAcc=zeros(opts.MaxIter,1);
    cntAcc=zeros(opts.MaxIter,1);
    fAcc(1) = arfitness;
    
    % Trace of best function values
    xBest=zeros(opts.MaxIter/opts.SavingModulo,N);
    xBest(1,:)=xstart;
    fBest=zeros(opts.MaxIter/opts.SavingModulo,1);
    fBest(1) = c_T;
    
    % Vector of step lengths
    rVec=zeros(opts.MaxIter/opts.SavingModulo,1);
    rVec(1)=r;
    
    % Vector of mu
    muVec = zeros(opts.MaxIter/opts.SavingModulo,N);
    muVec(1,:)=mu;
    
    % Vector of acceptance thresholds
    c_TVec=zeros(opts.MaxIter/opts.SavingModulo,1);
    c_TVec(1)=c_T;
    
    % Vector of empirical acceptance probability
    P_empVec = zeros(opts.MaxIter/opts.SavingModulo,1);
    P_empVec(1)=opts.valP;
    
    % Cell array of Q matrices
    if opts.bSaveCov==1
        QCell = cell(opts.MaxIter/opts.SavingModulo,1);
        QCell{1} = Q;
    end
    
    % Vector of evaluation indices
    countVec = zeros(opts.MaxIter/opts.SavingModulo,1);
    countVec(1) = 1;
    
    % Cell array with stop flags
    stopFlagCell = cell(1,1);
    
end

% Final initialization
saveInd = 2;
stopFlag='';

% -------------------- Generation Loop --------------------------------

% Reserve the last evaluation for the fitness of mu
while counteval < (opts.MaxIter-1)
    
    arfitness_old=arfitness;
    counteval = counteval+1;
    
    if strcmp(opts.Display,'on') && (mod(counteval,opts.VerboseModulo)==0)
        figure(23)
        error_ellipse(C(1:2,1:2),mu(1:2,1),'style','k')
        grid on
    end
    
    % Generate a N(0,1) variable
    arz = randn(N,1);  % array of normally distributed mutation vectors
    
    % Sample from N(mu,C) distribution
    arx_old = arx;
    arx = mu + r * (Q * arz);
    
    % Handle boundaries
    if opts.BoundActive
        arxNoBound=arx;
        arx=max([arx,opts.LBounds],[],2);
        arx=min([arx,opts.UBounds],[],2);
        if opts.BoundPenalty
            penalty=1e4*norm(arx-arxNoBound).^2;
        else
            penalty=0;
        end
    end
    
    
    % Plug-in local optimizer
    % arx=fminunc(fitfun,arx);
    
    % Objective function call
    if isempty(opts.funArgs)
        arfitness=feval(fitfun, arx)+penalty;
        % Just for recording
        % muFitness=feval(fitfun, mu);
    else
        arfitness=feval(fitfun, arx,opts.funArgs)+penalty;
        % Just for recording
        % muFitness=feval(fitfun, mu,opts.funArgs);
    end
    
    % Save best fitness values per restart
    if arfitness<currfBest
        currfBest = arfitness;
        currxBest = arx';
    end
    
    % Save best ever fitness values across restarts
    if arfitness<fmin
        fmin = arfitness;
        xmin = arx;
    end
    
    if strcmp(opts.Plotting,'on') && (mod(counteval,opts.VerboseModulo)==0) && strcmp(opts.bSaving,'on')
        
        if (opts.StopFitness  > -Inf)
            biasVal = -opts.StopFitness;
        else
            biasVal = 0;
        end
        
        currInd = counteval-1;
        currIndices = lastInd:opts.SavingModulo:currInd;
        currSaves = lastSave:saveInd-1;
        
        [eigVecs,eigVals] = eig(C./r^2); % eigen decomposition for plotting
        eigVals=diag(eigVals);
        
        try
            figure(42)
            subplot(2,2,1)
            plot([currIndices],muVec(currSaves,:)')
            title(['Current best: ',num2str(currfBest)])
            grid on
            hold on
            
            subplot(2,2,2)
            semilogy([currIndices],rVec(currSaves,:)')
            title(['Current step size: ',num2str(r)])
            grid on
            hold on
            
            subplot(2,2,3)
            semilogy([currIndices],biasVal+fBest([currSaves]))
            hold on
            grid on
            semilogy([currIndices],biasVal+c_TVec(currSaves,:)','r')
            title(['Current constraint: ',num2str(c_T)])
            
            %         subplot(2,2,4)
            %         plot([currIndices],P_empVec([currSaves]))
            %         grid on
            %         hold on
            %         title(['Current acceptance prob: ',num2str(P_emp)])
            
            subplot(2,2,4)
            semilogy([lastInd,currInd],[1./lastEigVals,1./eigVals],'-')
            grid on
            hold on
            title(['Condition of C: ',num2str(cond(C./r^2))])
        end
        lastEigVals = eigVals;
        
        lastInd = currInd;
        lastSave = saveInd-1;
        
        drawnow;
    end
    
    if (opts.Optim==1)
        % Check whether fitness/rank-based selection is used
        if opts.ThreshRank
            threshold = sortedHist(floor(opts.valP*histLen));
        else
            threshold = c_T;
        end
        
        % Check whether to adapt or not
        bAdapt = (arfitness<threshold);
        
    else
        % Boltzmann criterion for Adaptive MH
        accProb = min(1,arfitness/arfitness_old);
        % accProb = min(1,exp(arfitness_old-arfitness));
        bAdapt = (rand<accProb);
    end
    
    
    
    % Adapt moments upon acceptance
    if bAdapt
        
        % Count accepted solutions
        numAcc = numAcc + 1;
        rNumAcc = rNumAcc + 1;
        
        if strcmp(opts.bSaving,'on')
            cntAcc(numAcc,1)=counteval;
            xAcc(numAcc,:)=arx';
            fAcc(numAcc,1)=arfitness;
        end
        
        % Expand r by factor ss
        r = r*ss;
        
        % Adapt c_T
        if (opts.Optim==1)
            
            % Update history of accepted f and x
            histFAcc(2:end)=histFAcc(1:end-1);
            histFAcc(1)=arfitness;
            sortedHist = sort(histFAcc,'ascend');
            
            histXAcc(2:end,:)=histXAcc(1:end-1,:);
            histXAcc(1,:)=arx';
            
            if opts.ThreshRank
                
                % threshold for last accepted sample
                c_T = threshold;
                
            else
                
                % Standard exponential threshold decrease
                c_T = (1-1/N_T)*c_T + arfitness/N_T;
                
                % Linear threshold decrease
                % c_T = (1-counteval/opts.MaxIter)*300;
            end
        end
        
        % C_direct = (1-1/N_C)*C + ((mu-arx)*(mu-arx)')/N_C
        
        % Adapt mu
        mu_old=mu;
        mu = (1-1/N_mu)*mu + arx/N_mu;
        
        deltaC =  (1-1/N_C)*eye(N,N) + arz*arz'./N_C;
        % set deltaC to identity if covariance adaptation is not
        % required
        % deltaC =  eye(N,N);
        
        % Adapt Q
        
        deltaC=triu(deltaC)+triu(deltaC,1)'; % enforce symmetry
        [B,deltaQ] = eig(deltaC);       % eigen decomposition, B==normalized eigenvectors
        deltaQ = B*diag(sqrt(diag(deltaQ)))*B'; % deltaQ contains standard deviations now
        
        if any(imag(deltaQ(:)))
            deltaQ
            error('deltaQ is imaginary')
        end
        
        Q = Q*deltaQ;
        
        if any(imag(Q(:)))
            Q
            error('Q is imaginary')
        end
        
        % M_scaled = (r*Q)*(r*Q)'
        % M_before = (r/ss*Q)*(r/ss*Q)'
        
        % Normalize volume of Q to 1
        detQ=det(Q);
        Q = Q./((detQ).^(1/N));
        
        % C_after = (r/ss*Q)*(r/ss*Q)'
        
        % update C
        % C =  triu(C)+triu(C,1)'; % enforce symmetry
        
        % update Q
        % [B,Q] = eig(C);       % eigen decomposition, B==normalized eigenvectors
        % Q = B*diag(sqrt(diag(Q)))*B'; % D contains standard deviations now
        
        % detQ=det(Q);
        % Q = Q./((det(Q)).^(1/N));
        
    else
        
        % Contract r by factor sf
        r=r*sf;
        if (opts.Optim == 0)
            arfitness = arfitness_old;
            arx = arx_old;
        end
    end
    
    
    % Break, if fitness is good enough
    
    if arfitness <= opts.StopFitness
        fmin=arfitness;
        xmin=arx;
        stopFlag='fitness';
        break;
    end
    
    
    % Restart when any of the fitness, x or r tolerance are met
    
    %%%%% TolFun %%%%%
    currTolFun = abs(max(histFAcc)-min(histFAcc));
    if (currTolFun<opts.TolFun) && rNumAcc>histLen
        stopFlag='TolFun';
    end
    %     if (abs(arfitness_old-arfitness)<opts.TolFun)
    %         stopFlag='TolFun';
    %     end
    
    %%%%% TolX %%%%%
    currTolX = norm(histXAcc(1,:)-histXAcc(end,:));
    if (currTolX<opts.TolX)&& rNumAcc>histLen
        stopFlag='TolX';
    end
    
    %%%%% TolR %%%%%
    if (r<opts.TolR)
        stopFlag='TolR';
    end
    
    %%%%% TolCon %%%%%
    currTolCon = abs(c_T-currfBest);
    if currTolCon<opts.TolCon && rNumAcc>histLen
        stopFlag='TolCon';
    end
    
    % Upper bound for step size (hard-coded for now)
    r = min(r,rMax);
    
    % Restart when Restart flag is active AND stop flags are active
    if opts.bRestart && opts.Optim && ~(strcmp(stopFlag,''))
        
        % Show reason for restart
        disp([num2str(counteval),': Restart due to: ',stopFlag]);
        
        % Record stop flags
        if size(stopFlagCell,1)==1
            stopFlagCell{1,1}=stopFlag;
        else
            stopFlagCell{end+1,1}=stopFlag;
        end
        
        % Reset stopFlag
        stopFlag='';
        
        % Reset Gaussian parameters
        r = r_init;
        Q = eye(N,N);
        
        % Handle boundaries
        if opts.BoundActive==1
            mu = opts.LBounds+rand(N,1).*(opts.UBounds-opts.LBounds);
            arx = mu;
            
        else
            mu = xstart;
            arx = xstart;
        end
        
        % New start value
        if isempty(opts.funArgs)
            c_T=feval(fitfun, arx);
        else
            c_T=feval(fitfun, arx,opts.funArgs);
        end
        
        counteval = counteval+1;
        arfitness = c_T;
        rNumAcc = 0;
        
        % Restart with larger N_T
        N_T = inc_T*N_T;
        
        % Increaed history length for rank-based selection
        histFAcc = arfitness*ones(histLen,1);
        sortedHist = histFAcc;
        
        currfBest=arfitness;
        histXAcc = repmat(arx',histLen,1);
        
    end
    
    C =  (r*Q)*(r*Q)';
    
    P_emp = numAcc/counteval;
    
    
    % Save data
    if strcmp(opts.bSaving,'on') && (mod(counteval,opts.SavingModulo)==0)
        
        rVec(saveInd)=r;
        muVec(saveInd,:) = mu;
        c_TVec(saveInd)=c_T;
        P_empVec(saveInd)=P_emp;
        % Cell array of Q matrices
        
        if opts.bSaveCov==1
            
            QCell{saveInd}=Q;
        end
        %if arfitness<=threshold
        fRaw(saveInd)=arfitness;
        xRaw(saveInd,:)=arx;
        %else
        %    fRaw(saveInd)=fRaw(saveInd-1);
        %    xRaw(saveInd,:)=xRaw(saveInd-1,:);
        %end
        fBest(saveInd)=currfBest;
        xBest(saveInd,:)=currxBest;
        
        countVec(saveInd)=counteval;
        
        saveInd = saveInd+1;
        
    end
    
    if (mod(counteval,opts.VerboseModulo)==0)
        disp( '-------------------------------------------');
        disp([' Number of iterations: ',num2str(counteval)]);
        disp([' Current fitness:      ',num2str(arfitness)]);
        disp([' Constraint:           ',num2str(c_T)]);
        disp([' P_acc:                ',num2str(numAcc/counteval)]);
        disp([' Search radius:        ',num2str(r)]);
        disp([' TolCon:               ',num2str(currTolCon)])
        disp([' TolX:                 ',num2str(currTolX)])
    end
    
end % while, end generation loop

if (counteval==(opts.MaxIter-1))
    
    % Final objective function call is reserved for mu
    if isempty(opts.funArgs)
        muFitness=feval(fitfun, mu);
    else
        muFitness=feval(fitfun, mu,opts.funArgs);
    end
    
    if muFitness <= opts.StopFitness
        currfBest = muFitness;
        currxBest = mu;
        stopFlag='fitness';
        if muFitness <=fmin
            fmin = muFitness;
            xmin = mu;
        end
    else
        stopFlag='MaxIter';
    end
    counteval = counteval+1;
end

if strcmp(opts.bSaving,'on')
    
    % Save the last evaluation
    rVec(saveInd)=r;
    muVec(saveInd,:) = mu;
    xRaw(saveInd,:) = arx;
    fRaw(saveInd,:) = arfitness;
    c_TVec(saveInd)=c_T;
    P_empVec(saveInd)=P_emp;
    
    if opts.bSaveCov==1
        QCell{saveInd}=Q;
    end
    
    fBest(saveInd)=currfBest;
    xBest(saveInd,:)=currxBest;
    countVec(saveInd)=counteval;
    
    out.xRaw = xRaw(1:saveInd,:);
    out.fRaw = fRaw(1:saveInd,:);
    out.cntAcc=cntAcc(1:numAcc,:);
    out.fAcc=fAcc(1:numAcc,:);
    out.xAcc=xAcc(1:numAcc,:);
    out.xBest = xBest(1:saveInd,:);
    out.fBest = fBest(1:saveInd,:);
    out.P_empVec = P_empVec(1:saveInd,:);
    out.rVec = rVec(1:saveInd,:);
    out.muVec = muVec(1:saveInd,:);
    out.c_TVec = c_TVec(1:saveInd,:);
    
    out.bestever.x = xmin;
    out.bestever.f = fmin;
    
    if opts.bSaveCov==1
        out.QCell = QCell(1:saveInd,:);
    end
    out.countVec = countVec(1:saveInd);
    
    % Record all stop flags
    if isempty(stopFlagCell{1,1})
        stopFlagCell{1,1}=stopFlag;
    else
        stopFlagCell{end+1,1}=stopFlag;
    end
    out.stopFlagCell = stopFlagCell;
    
end

out.stopFlag = stopFlag;
out.P_emp = numAcc/counteval;
out.opts = opts;

% -------------------- Ending Message ----------------------------

disp(['          Acceptance probability: ' num2str(numAcc/counteval)]);
if (opts.Optim == 1)
    disp([num2str(counteval) ' Current fitness: ' num2str(arfitness)]);
    disp(['          Best ever fitness: ' num2str(fmin)]);
end
disp(['          Stop flag: ' stopFlag]);

if strcmp(opts.Display,'on')
    error_ellipse('C',C(1:2,1:2),'mu',mu(1:2),'style','b')
end

% ------------------- Test functions -----------------------------
% PARTLY TAKEN FROM CMA-ES
% ----------------------------------------------------------------

function f=fsphere(x)
f=sum(x.^2);

function f=flogsphere(x)
f=log(1+sum(x.^2));

function f=fschwefel(x)
f = sum(cumsum(x).^2);

function f=fcigar(x)
f = x(1)^2 + 1e6*sum(x(2:end).^2);

function f=fcigtab(x)
f = x(1)^2 + 1e8*x(end)^2 + 1e4*sum(x(2:(end-1)).^2);

function f=ftablet(x)
f = 1e6*x(1)^2 + sum(x(2:end).^2);

function f=felli(x)
N = size(x,1); if N < 2 error('dimension must be greater one'); end
f=1e6.^((0:N-1)/(N-1)) * x.^2;

function f=felli100(x)
N = size(x,1); if N < 2 error('dimension must be greater one'); end
f=1e4.^((0:N-1)/(N-1)) * x.^2;

function f=fplane(x)
f=x(1);

function f=fBaseball(x)
f=sqrt(sum(x.^2));

function f=frosen(x)
if size(x,1) < 2 error('dimension must be greater one'); end
f = 100*sum((x(1:end-1).^2 - x(2:end)).^2) + sum((x(1:end-1)-1).^2);


function f=frastrigin(x)
N = size(x,1); if N < 2 error('dimension must be greater one'); end
amp=10;
f = amp * N + sum(x .^2 - amp * cos(2 * pi .* x),1);

function f=flograstrigin(x)
% Check for scale-invariance
N = size(x,1); if N < 2 error('dimension must be greater one'); end
amp=10;
f = amp * N + sum(x .^2 - amp * cos(2 * pi .* x),1);
f=log(1+f);

function f=fdsphere(x,params)
% Lunacek et al. 2008
s=params(1);
d=params(2);
N=length(x);
mu1=2.5*ones(N,1);
mu2=-sqrt((mu1.^2-d)./s);
f=min(sum((x-mu1).^2),d*N+s*sum((x-mu2).^2));

function f=fdrastrigin(x,params)
% Lunacek et al. 2008
N = length(x);
mu1 = 2.5*ones(N,1);

if length(params)==2
    s=params(1);
    d=params(2);
    mu2=-sqrt((mu1.^2-d)./s);
else
    s=params(1);
    mu2 = -2.5*ones(N,1);
    d = mu1(1).^2-s*mu2(1).^2;
end

f1=min(sum((x-mu1).^2),d*N+s*sum((x-mu2).^2));
x_mu=(x-mu1);
f2= 10*sum(1-cos(2 * pi .*x_mu),1);
f=f1+f2;

function f=fnoisysphere(x,sigma_eps)
% Additive noise
randNum = sigma_eps*randn;
f=norm(x)^2+randNum;

% Multiplicative noise
% randNum = -1.5+3*rand;
% f=norm(x)^2*(1+randNum);


function f=frand(x)
f=rand;


%%%%%%%%% Examples from Kjellstrom 1981-1999 %%%%%%%%%

function f=fkjellstrom1(x)
N = size(x,1); if N < 2 error('dimension must be greater one'); end
f = 20.0 - sum((x+0.5).^2) + exp(25*(sum(x.^2)-8));

function f=fkjellstrom2(x)
N = size(x,1); if N < 2 error('dimension must be greater one'); end
h = 0.01*(cos(x+1.982) + cos(2*x+5.720) + cos(3*x + 1.621) + cos(4*x + 0.823) + cos(5*x + 3.22));
f = prod(1+h(:));


%%%%%%%%% Examples for test target distributions %%%%%%%%%

function f=fHaario(x,type)
% from Haario et al. Comp Stat. 1999
global HaarioMat;
N = size(x,1);
C=eye(N,N);
switch type
    case 1
        C(1,1)=100;
    case 2
        C(1,1)=100;
        % TODO rotation
    case 3
        C(1,1)=100;
        b=0.03;
        x(2,1)=x(2,1)+b*x(1,1)^2-100*b;
    case 4
        C(1,1)=100;
        b=0.1;
        x(2,1)=x(2,1)+b*x(1,1)^2-100*b;
    otherwise
end
f = -mvnpdf(x,zeros(N,1),C);

function f=fLiangEx1(x)
% taken from Liang:2001
global LiangExMat;
sigma=0.1;
w=0.05;
if length(x)~=2
    error('Wrong dimension')
else
    f=0;
    for i=1:2:39
        currCov = (x-LiangExMat(i:i+1))'*(x-LiangExMat(i:i+1));
        f=f+w*exp(-1/(2*sigma^2)*currCov);
    end
    f=-1/sqrt(2*pi*sigma)*f;
end

function f=fLiangEx2(x)
% taken from Liang:2001
N = length(x);
f = -1/3*mvnpdf(x,zeros(N,1),eye(N,N));
f = f - 2/3*mvnpdf(x,5*ones(N,1),eye(N,N));

% ---------------------------------------------------------------
% FUNCTIONS BELOW ARE TAKEN FROM Niko Hansen's CMA-ES code
% ---------------------------------------------------------------
function opts=getoptions(inopts, defopts)
% OPTS = GETOPTIONS(INOPTS, DEFOPTS) handles an arbitrary number of
% optional arguments to a function. The given arguments are collected
% in the struct INOPTS.  GETOPTIONS matches INOPTS with a default
% options struct DEFOPTS and returns the merge OPTS.  Empty or missing
% fields in INOPTS invoke the default value.  Fieldnames in INOPTS can
% be abbreviated.
%
% The returned struct OPTS is first assigned to DEFOPTS. Then any
% field value in OPTS is replaced by the respective field value of
% INOPTS if (1) the field unambiguously (case-insensitive) matches
% with the fieldname in INOPTS (cut down to the length of the INOPTS
% fieldname) and (2) the field is not empty.
%
% Example:
%   In the source-code of the function that needs optional
%   arguments, the last argument is the struct of optional
%   arguments:
%
%   function results = myfunction(mandatory_arg, inopts)
%     % Define four default options
%     defopts.PopulationSize = 200;
%     defopts.ParentNumber = 50;
%     defopts.MaxIterations = 1e6;
%     defopts.MaxSigma = 1;
%
%     % merge default options with input options
%     opts = getoptions(inopts, defopts);
%
%     % Thats it! From now on the values in opts can be used
%     for i = 1:opts.PopulationSize
%       % do whatever
%       if sigma > opts.MaxSigma
%         % do whatever
%       end
%     end
%
%   For calling the function myfunction with default options:
%   myfunction(argument1, []);
%   For calling the function myfunction with modified options:
%   opt.pop = 100; % redefine PopulationSize option
%   opt.PAR = 10;  % redefine ParentNumber option
%   opt.maxiter = 2; % opt.max=2 is ambiguous and would result in an error
%   myfunction(argument1, opt);

%
% 04/07/19: Entries can be structs itself leading to a recursive
%           call to getoptions.
%

if nargin < 2 || isempty(defopts) % no default options available
    opts=inopts;
    return;
elseif isempty(inopts) % empty inopts invoke default options
    opts = defopts;
    return;
elseif ~isstruct(defopts) % handle a single option value
    if isempty(inopts)
        opts = defopts;
    elseif ~isstruct(inopts)
        opts = inopts;
    else
        error('Input options are a struct, while default options are not');
    end
    return;
elseif ~isstruct(inopts) % no valid input options
    error('The options need to be a struct or empty');
end

opts = defopts; % start from defopts
% if necessary overwrite opts fields by inopts values
defnames = fieldnames(defopts);
idxmatched = []; % indices of defopts that already matched
for name = fieldnames(inopts)'
    name = name{1}; % name of i-th inopts-field
    if isoctave
        for i = 1:size(defnames, 1)
            idx(i) = strncmp(lower(defnames(i)), lower(name), length(name));
        end
    else
        idx = strncmp(lower(defnames), lower(name), length(name));
    end
    if sum(idx) > 1
        error(['option "' name '" is not an unambigous abbreviation. ' ...
            'Use opts=RMFIELD(opts, ''' name, ...
            ''') to remove the field from the struct.']);
    end
    if sum(idx) == 1
        defname  = defnames{find(idx)};
        if ismember(find(idx), idxmatched)
            error(['input options match more than ones with "' ...
                defname '". ' ...
                'Use opts=RMFIELD(opts, ''' name, ...
                ''') to remove the field from the struct.']);
        end
        idxmatched = [idxmatched find(idx)];
        val = getfield(inopts, name);
        % next line can replace previous line from MATLAB version 6.5.0 on and in octave
        % val = inopts.(name);
        if isstruct(val) % valid syntax only from version 6.5.0
            opts = setfield(opts, defname, ...
                getoptions(val, getfield(defopts, defname)));
        elseif isstruct(getfield(defopts, defname))
            % next three lines can replace previous three lines from MATLAB
            % version 6.5.0 on
            %   opts.(defname) = ...
            %      getoptions(val, defopts.(defname));
            % elseif isstruct(defopts.(defname))
            warning(['option "' name '" disregarded (must be struct)']);
        elseif ~isempty(val) % empty value: do nothing, i.e. stick to default
            opts = setfield(opts, defnames{find(idx)}, val);
            % next line can replace previous line from MATLAB version 6.5.0 on
            % opts.(defname) = inopts.(name);
        end
    else
        warning(['option "' name '" disregarded (unknown field name)']);
    end
end


% ---------------------------------------------------------------
% ---------------------------------------------------------------
function res = isoctave
% any hack to find out whether we are running octave
s = version;
res = 0;
if exist('fflush', 'builtin') && eval(s(1)) < 7
    res = 1;
end