• Main Page
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

BCIntegrate.cxx

Go to the documentation of this file.
00001 /*
00002  * Copyright (C) 2008-2010, Daniel Kollar and Kevin Kroeninger.
00003  * All rights reserved.
00004  *
00005  * For the licensing terms see doc/COPYING.
00006  */
00007 
00008 // ---------------------------------------------------------
00009 
00010 #include "config.h"
00011 
00012 #include "BAT/BCIntegrate.h"
00013 #include "BAT/BCLog.h"
00014 #include "BAT/BCMath.h"
00015 
00016 #include <TH1D.h>
00017 #include <TH2D.h>
00018 #include <TMinuit.h>
00019 #include <TRandom3.h>
00020 #include <TTree.h>
00021 
00022 #ifdef HAVE_CUBA_H
00023 #include "cuba.h"
00024 #endif
00025 
00026 // ---------------------------------------------------------
00027 
00028 class BCIntegrate * global_this;
00029 
00030 // *********************************************
00031 BCIntegrate::BCIntegrate() : BCEngineMCMC()
00032 {
00033    fNvar=0;
00034    fNiterPerDimension = 100;
00035    fNSamplesPer2DBin = 100;
00036 
00037    fNIterationsMax    = 1000000;
00038    fRelativePrecision = 1e-3;
00039 
00040 #ifdef HAVE_CUBA_H
00041    fIntegrationMethod = BCIntegrate::kIntCuba;
00042 #else
00043    fIntegrationMethod = BCIntegrate::kIntMonteCarlo;
00044 #endif
00045    fMarginalizationMethod = BCIntegrate::kMargMetropolis;
00046    fOptimizationMethod = BCIntegrate::kOptMinuit;
00047    fOptimizationMethodMode = BCIntegrate::kOptMinuit;
00048 
00049    fNbins=100;
00050 
00051    fDataPointLowerBoundaries = 0;
00052    fDataPointUpperBoundaries = 0;
00053 
00054    fFillErrorBand = false;
00055 
00056    fFitFunctionIndexX = -1;
00057    fFitFunctionIndexY = -1;
00058 
00059    fErrorBandNbinsX = 100;
00060    fErrorBandNbinsY = 500;
00061 
00062    fErrorBandContinuous = true;
00063 
00064    fMinuit = 0;
00065    fMinuitArglist[0] = 20000;
00066    fMinuitArglist[1] = 0.01;
00067    fFlagIgnorePrevOptimization = false;
00068 
00069    fFlagWriteMarkovChain = false;
00070    fMarkovChainTree = 0;
00071 
00072    fMarkovChainStepSize = 0.1;
00073 
00074    fMarkovChainAutoN = true;
00075 
00076    fSAT0 = 100.0;
00077    fSATmin = 0.1;
00078    fSASchedule = BCIntegrate::kSACauchy;
00079    fFlagWriteSAToFile = false;
00080 }
00081 
00082 // *********************************************
00083 BCIntegrate::BCIntegrate(BCParameterSet * par) : BCEngineMCMC(1)
00084 {
00085    fNvar=0;
00086    fNiterPerDimension = 100;
00087    fNSamplesPer2DBin = 100;
00088 
00089    fNbins=100;
00090 
00091    SetParameters(par);
00092 
00093    fDataPointLowerBoundaries = 0;
00094    fDataPointUpperBoundaries = 0;
00095 
00096    fFillErrorBand = false;
00097 
00098    fFitFunctionIndexX = -1;
00099    fFitFunctionIndexY = -1;
00100 
00101    fErrorBandNbinsX = 100;
00102    fErrorBandNbinsY = 200;
00103 
00104    fErrorBandContinuous = true;
00105 
00106    fMinuit = 0;
00107    fMinuitArglist[0] = 20000;
00108    fMinuitArglist[1] = 0.01;
00109    fFlagIgnorePrevOptimization = false;
00110 
00111    fFlagWriteMarkovChain = false;
00112    fMarkovChainTree = 0;
00113 
00114    fMarkovChainStepSize = 0.1;
00115 
00116    fMarkovChainAutoN = true;
00117 
00118    fSAT0 = 100.0;
00119    fSATmin = 0.1;
00120    fTreeSA = false;
00121 }
00122 
00123 // *********************************************
00124 BCIntegrate::~BCIntegrate()
00125 {
00126    DeleteVarList();
00127 
00128    fx=0;
00129 
00130    if (fMinuit)
00131       delete fMinuit;
00132 
00133    int n1 = fHProb1D.size();
00134    if(n1>0) {
00135       for (int i=0;i<n1;i++)
00136          delete fHProb1D.at(i);
00137    }
00138 
00139    int n2 = fHProb2D.size();
00140    if(n2>0) {
00141       for (int i=0;i<n2;i++)
00142          delete fHProb2D.at(i);
00143    }
00144 }
00145 
00146 // *********************************************
00147 void BCIntegrate::SetParameters(BCParameterSet * par)
00148 {
00149    DeleteVarList();
00150 
00151    fx = par;
00152    fNvar = fx->size();
00153    fMin = new double[fNvar];
00154    fMax = new double[fNvar];
00155    fVarlist = new int[fNvar];
00156 
00157    ResetVarlist(1);
00158 
00159    for(int i=0;i<fNvar;i++) {
00160       fMin[i]=(fx->at(i))->GetLowerLimit();
00161       fMax[i]=(fx->at(i))->GetUpperLimit();
00162    }
00163 
00164    fXmetro0.clear();
00165    for(int i=0;i<fNvar;i++)
00166       fXmetro0.push_back((fMin[i]+fMax[i])/2.0);
00167 
00168    fXmetro1.clear();
00169    fXmetro1.assign(fNvar,0.);
00170 
00171    fMCMCBoundaryMin.clear();
00172    fMCMCBoundaryMax.clear();
00173 
00174    for(int i=0;i<fNvar;i++) {
00175       fMCMCBoundaryMin.push_back(fMin[i]);
00176       fMCMCBoundaryMax.push_back(fMax[i]);
00177       fMCMCFlagsFillHistograms.push_back(true);
00178    }
00179 
00180    for (int i = int(fMCMCH1NBins.size()); i<fNvar; ++i)
00181       fMCMCH1NBins.push_back(100);
00182 
00183    fMCMCNParameters = fNvar;
00184 }
00185 
00186 // *********************************************
00187 void BCIntegrate::SetMarkovChainInitialPosition(std::vector<double> position)
00188 {
00189    int n = position.size();
00190 
00191    fXmetro0.clear();
00192 
00193    for (int i = 0; i < n; ++i)
00194       fXmetro0.push_back(position.at(i));
00195 }
00196 
00197 // *********************************************
00198 void BCIntegrate::DeleteVarList()
00199 {
00200    if(fNvar) {
00201       delete[] fVarlist;
00202       fVarlist=0;
00203 
00204       delete[] fMin;
00205       fMin=0;
00206 
00207       delete[] fMax;
00208       fMax=0;
00209 
00210       fx=0;
00211       fNvar=0;
00212    }
00213 }
00214 
00215 // *********************************************
00216 void BCIntegrate::SetNbins(int nbins, int index)
00217 {
00218    if (fNvar == 0)
00219       return;
00220 
00221    // check if index is in range
00222    if (index >= fNvar) {
00223       BCLog::OutWarning("BCIntegrate::SetNbins : Index out of range.");
00224       return;
00225    }
00226    // set for all parameters at once
00227    else if (index < 0) {
00228       for (int i = 0; i < fNvar; ++i)
00229          SetNbins(nbins, i);
00230       return;
00231    }
00232 
00233    // sanity check for nbins
00234    if (nbins <= 0)
00235       nbins = 10;
00236 
00237    fMCMCH1NBins[index] = nbins;
00238 
00239    return;
00240 
00241 //    if(n<2) {
00242 //       BCLog::OutWarning("BCIntegrate::SetNbins. Number of bins less than 2 makes no sense. Setting to 2.");
00243 //       n=2;
00244 //    }
00245 //    MCMCSetH1NBins(n, -1);
00246 
00247    //   fNbins=n;
00248 
00249    //   fMCMCH1NBins = n;
00250    //   fMCMCH2NBinsX = n;
00251    //   fMCMCH2NBinsY = n;
00252 }
00253 
00254 // *********************************************
00255 // void BCIntegrate::SetNbinsX(int n)
00256 // {
00257 //    if(n<2) {
00258 //       BCLog::OutWarning("BCIntegrate::SetNbins. Number of bins less than 2 makes no sense. Setting to 2.");
00259 //       n=2;
00260 //    }
00261 //    fMCMCH2NBinsX = n;
00262 // }
00263 
00264 // *********************************************
00265 // void BCIntegrate::SetNbinsY(int n)
00266 // {
00267 //    if(n<2) {
00268 //       BCLog::OutWarning("BCIntegrate::SetNbins. Number of bins less than 2 makes no sense. Setting to 2.");
00269 //       n=2;
00270 //    }
00271 //    fNbins=n;
00272 
00273 //    fMCMCH2NBinsY = n;
00274 // }
00275 
00276 // *********************************************
00277 void BCIntegrate::SetVarList(int * varlist)
00278 {
00279    for(int i=0;i<fNvar;i++)
00280       fVarlist[i]=varlist[i];
00281 }
00282 
00283 // *********************************************
00284 void BCIntegrate::ResetVarlist(int v)
00285 {
00286    for(int i=0;i<fNvar;i++)
00287       fVarlist[i]=v;
00288 }
00289 
00290 // *********************************************
00291 double BCIntegrate::Eval(std::vector <double> x)
00292 {
00293    BCLog::OutWarning( "BCIntegrate::Eval. No function. Function needs to be overloaded.");
00294    return 0;
00295 }
00296 
00297 // *********************************************
00298 double BCIntegrate::LogEval(std::vector <double> x)
00299 {
00300    // this method should better also be overloaded
00301    return log(Eval(x));
00302 }
00303 
00304 // *********************************************
00305 double BCIntegrate::EvalSampling(std::vector <double> x)
00306 {
00307    BCLog::OutWarning( "BCIntegrate::EvalSampling. No function. Function needs to be overloaded.");
00308    return 0;
00309 }
00310 
00311 // *********************************************
00312 double BCIntegrate::LogEvalSampling(std::vector <double> x)
00313 {
00314    return log(EvalSampling(x));
00315 }
00316 
00317 // *********************************************
00318 double BCIntegrate::EvalPrint(std::vector <double> x)
00319 {
00320    double val=Eval(x);
00321    BCLog::OutDebug(Form("BCIntegrate::EvalPrint. Value: %d.", val));
00322 
00323    return val;
00324 }
00325 
00326 // *********************************************
00327 void BCIntegrate::SetIntegrationMethod(BCIntegrate::BCIntegrationMethod method)
00328 {
00329 #ifdef HAVE_CUBA_H
00330    fIntegrationMethod = method;
00331 #else
00332    BCLog::OutWarning("!!! This version of BAT is compiled without Cuba.");
00333    BCLog::OutWarning("    Monte Carlo Sampled Mean integration method will be used.");
00334    BCLog::OutWarning("    To be able to use Cuba integration, install Cuba and recompile BAT.");
00335 #endif
00336 }
00337 
00338 // *********************************************
00339 int BCIntegrate::IntegrateResetResults()
00340 {
00341    fBestFitParameters.clear();
00342    fBestFitParameterErrors.clear();
00343    fBestFitParametersMarginalized.clear();
00344 
00345    // no errors
00346    return 1;
00347 }
00348 
00349 // *********************************************
00350 double BCIntegrate::Integrate()
00351 {
00352    std::vector <double> parameter;
00353    parameter.assign(fNvar, 0.0);
00354 
00355    switch(fIntegrationMethod)
00356    {
00357       case BCIntegrate::kIntMonteCarlo:
00358          return IntegralMC(parameter);
00359 
00360       case BCIntegrate::kIntMetropolis:
00361          return IntegralMetro(parameter);
00362 
00363       case BCIntegrate::kIntImportance:
00364          return IntegralImportance(parameter);
00365 
00366       case BCIntegrate::kIntCuba:
00367 #ifdef HAVE_CUBA_H
00368          return CubaIntegrate();
00369 #else
00370          BCLog::OutError("!!! This version of BAT is compiled without Cuba.");
00371          BCLog::OutError("    Use other integration methods or install Cuba and recompile BAT.");
00372 #endif
00373    }
00374 
00375    BCLog::OutError(
00376          Form("BCIntegrate::Integrate : Invalid integration method: %d", fIntegrationMethod));
00377 
00378    return 0;
00379 }
00380 
00381 // *********************************************
00382 double BCIntegrate::IntegralMC(std::vector <double> x, int * varlist)
00383 {
00384    SetVarList(varlist);
00385    return IntegralMC(x);
00386 }
00387 
00388 // *********************************************
00389 double BCIntegrate::IntegralMC(std::vector <double> x)
00390 {
00391    // count the variables to integrate over
00392    int NvarNow=0;
00393 
00394    for(int i = 0; i < fNvar; i++)
00395       if(fVarlist[i])
00396          NvarNow++;
00397 
00398    // print to log
00399    BCLog::OutDebug(Form(" --> Running MC integation over %i dimensions.", NvarNow));
00400    BCLog::OutDebug(Form(" --> Maximum number of iterations: %i", GetNIterationsMax()));
00401    BCLog::OutDebug(Form(" --> Aimed relative precision:     %e", GetRelativePrecision()));
00402 
00403    // calculate D (the integrated volume)
00404    double D = 1.0;
00405    for(int i = 0; i < fNvar; i++)
00406       if (fVarlist[i])
00407          D = D * (fMax[i] - fMin[i]);
00408 
00409    // reset variables
00410    double pmax = 0.0;
00411    double sumW  = 0.0;
00412    double sumW2 = 0.0;
00413    double integral = 0.0;
00414    double variance = 0.0;
00415    double relprecision = 1.0;
00416 
00417    std::vector <double> randx;
00418    randx.assign(fNvar, 0.0);
00419 
00420    // reset number of iterations
00421    fNIterations = 0;
00422 
00423    // iterate while precision is not reached and the number of iterations is lower than maximum number of iterations
00424    while ((fRelativePrecision < relprecision && fNIterationsMax > fNIterations) || fNIterations < 10) {
00425       // increase number of iterations
00426       fNIterations++;
00427 
00428       // get random numbers
00429       GetRandomVector(randx);
00430 
00431       // scale random numbers
00432       for(int i = 0; i < fNvar; i++) {
00433          if(fVarlist[i])
00434             randx[i]=fMin[i]+randx[i]*(fMax[i]-fMin[i]);
00435          else
00436             randx[i]=x[i];
00437       }
00438 
00439       // evaluate function at sampled point
00440       double value = Eval(randx);
00441 
00442       // add value to sum and sum of squares
00443       sumW  += value;
00444       sumW2 += value * value;
00445 
00446       // search for maximum probability
00447       if (value > pmax) {
00448          // set new maximum value
00449          pmax = value;
00450 
00451          // delete old best fit parameter values
00452          fBestFitParameters.clear();
00453 
00454          // write best fit parameters
00455          for(int i = 0; i < fNvar; i++)
00456             fBestFitParameters.push_back(randx.at(i));
00457       }
00458 
00459       // calculate integral and variance
00460       integral = D * sumW / fNIterations;
00461 
00462       if (fNIterations%10000 == 0) {
00463          variance = (1.0 / double(fNIterations)) * (D * D * sumW2 / double(fNIterations) - integral * integral);
00464          double error = sqrt(variance);
00465          relprecision = error / integral;
00466 
00467          BCLog::Out(BCLog::debug, BCLog::debug,
00468             Form("BCIntegrate::IntegralMC. Iteration %i, integral: %e +- %e.", fNIterations, integral, error));
00469       }
00470    }
00471 
00472    fError = variance / fNIterations;
00473 
00474    // print to log
00475    BCLog::OutDebug(Form(" --> Result of integration:        %e +- %e   in %i iterations.", integral, sqrt(variance), fNIterations));
00476    BCLog::OutDebug(Form(" --> Obtained relative precision:  %e. ", sqrt(variance) / integral));
00477 
00478    return integral;
00479 }
00480 
00481 
00482 // *********************************************
00483 double BCIntegrate::IntegralMetro(std::vector <double> x)
00484 {
00485    // print debug information
00486    BCLog::OutDebug(Form("BCIntegrate::IntegralMetro. Integate over %i dimensions.", fNvar));
00487 
00488    // get total number of iterations
00489    double Niter = pow(fNiterPerDimension, fNvar);
00490 
00491    // print if more than 100,000 iterations
00492    if(Niter>1e5)
00493       BCLog::OutDebug(Form(" --> Total number of iterations in Metropolis: %d.", Niter));
00494 
00495    // reset sum
00496    double sumI = 0;
00497 
00498    // prepare Metropolis algorithm
00499    std::vector <double> randx;
00500    randx.assign(fNvar,0.);
00501    InitMetro();
00502 
00503    // prepare maximum value
00504    double pmax = 0.0;
00505 
00506    // loop over iterations
00507    for(int i = 0; i <= Niter; i++) {
00508       // get random point from sampling distribution
00509       GetRandomPointSamplingMetro(randx);
00510 
00511       // calculate probability at random point
00512       double val_f = Eval(randx);
00513 
00514       // calculate sampling distributions at that point
00515       double val_g = EvalSampling(randx);
00516 
00517       // add ratio to sum
00518       if (val_g > 0)
00519          sumI += val_f / val_g;
00520 
00521       // search for maximum probability
00522       if (val_f > pmax) {
00523          // set new maximum value
00524          pmax = val_f;
00525 
00526          // delete old best fit parameter values
00527          fBestFitParameters.clear();
00528 
00529          // write best fit parameters
00530          for(int i = 0; i < fNvar; i++)
00531             fBestFitParameters.push_back(randx.at(i));
00532       }
00533 
00534       // write intermediate results
00535       if((i+1)%100000 == 0)
00536          BCLog::OutDebug(Form("BCIntegrate::IntegralMetro. Iteration %d, integral: %d.", i+1, sumI/(i+1)));
00537    }
00538 
00539    // calculate integral
00540    double result=sumI/Niter;
00541 
00542    // print debug information
00543    BCLog::OutDebug(Form(" --> Integral is %f in %i iterations.", result, Niter));
00544 
00545    return result;
00546 }
00547 
00548 // *********************************************
00549 double BCIntegrate::IntegralImportance(std::vector <double> x)
00550 {
00551    // print debug information
00552    BCLog::OutDebug(Form("BCIntegrate::IntegralImportance. Integate over %i dimensions.", fNvar));
00553 
00554    // get total number of iterations
00555    double Niter = pow(fNiterPerDimension, fNvar);
00556 
00557    // print if more than 100,000 iterations
00558    if(Niter>1e5)
00559       BCLog::OutDebug(Form("BCIntegrate::IntegralImportance. Total number of iterations: %d.", Niter));
00560 
00561    // reset sum
00562    double sumI = 0;
00563 
00564    std::vector <double> randx;
00565    randx.assign(fNvar,0.);
00566 
00567    // prepare maximum value
00568    double pmax = 0.0;
00569 
00570    // loop over iterations
00571    for(int i = 0; i <= Niter; i++) {
00572       // get random point from sampling distribution
00573       GetRandomPointImportance(randx);
00574 
00575       // calculate probability at random point
00576       double val_f = Eval(randx);
00577 
00578       // calculate sampling distributions at that point
00579       double val_g = EvalSampling(randx);
00580 
00581       // add ratio to sum
00582       if (val_g > 0)
00583          sumI += val_f / val_g;
00584 
00585       // search for maximum probability
00586       if (val_f > pmax) {
00587          // set new maximum value
00588          pmax = val_f;
00589 
00590          // delete old best fit parameter values
00591          fBestFitParameters.clear();
00592 
00593          // write best fit parameters
00594          for(int i = 0; i < fNvar; i++)
00595             fBestFitParameters.push_back(randx.at(i));
00596       }
00597 
00598       // write intermediate results
00599       if((i+1)%100000 == 0)
00600          BCLog::OutDebug(Form("BCIntegrate::IntegralImportance. Iteration %d, integral: %d.", i+1, sumI/(i+1)));
00601    }
00602 
00603    // calculate integral
00604    double result=sumI/Niter;
00605 
00606    // print debug information
00607    BCLog::OutDebug(Form("BCIntegrate::IntegralImportance. Integral %f in %i iterations.", result, Niter));
00608 
00609    return result;
00610 }
00611 
00612 // *********************************************
00613 TH1D* BCIntegrate::Marginalize(BCParameter * parameter)
00614 {
00615    BCLog::OutDebug(Form(" --> Marginalizing model wrt. parameter %s using method %d.", parameter->GetName().data(), fMarginalizationMethod));
00616 
00617    switch(fMarginalizationMethod)
00618    {
00619       case BCIntegrate::kMargMonteCarlo:
00620          return MarginalizeByIntegrate(parameter);
00621 
00622       case BCIntegrate::kMargMetropolis:
00623          return MarginalizeByMetro(parameter);
00624    }
00625 
00626    BCLog::OutWarning(Form("BCIntegrate::Marginalize. Invalid marginalization method: %d. Return 0.", fMarginalizationMethod));
00627 
00628    return 0;
00629 }
00630 
00631 // *********************************************
00632 TH2D * BCIntegrate::Marginalize(BCParameter * parameter1, BCParameter * parameter2)
00633 {
00634    switch(fMarginalizationMethod)
00635    {
00636       case BCIntegrate::kMargMonteCarlo:
00637          return MarginalizeByIntegrate(parameter1,parameter2);
00638 
00639       case BCIntegrate::kMargMetropolis:
00640          return MarginalizeByMetro(parameter1,parameter2);
00641    }
00642 
00643    BCLog::OutWarning(Form("BCIntegrate::Marginalize. Invalid marginalization method: %d. Return 0.", fMarginalizationMethod));
00644 
00645    return 0;
00646 }
00647 
00648 // *********************************************
00649 TH1D* BCIntegrate::MarginalizeByIntegrate(BCParameter * parameter)
00650 {
00651    // set parameter to marginalize
00652    ResetVarlist(1);
00653    int index = parameter->GetIndex();
00654    UnsetVar(index);
00655 
00656    // define histogram
00657    double hmin = parameter->GetLowerLimit();
00658    double hmax = parameter->GetUpperLimit();
00659    double hdx  = (hmax - hmin) / double(fNbins);
00660    TH1D * hist = new TH1D("hist","", fNbins, hmin, hmax);
00661 
00662    // integrate
00663    std::vector <double> randx;
00664    randx.assign(fNvar, 0.);
00665 
00666    for(int i=0;i<=fNbins;i++) {
00667       randx[index] = hmin + (double)i * hdx;
00668 
00669       double val = IntegralMC(randx);
00670       hist->Fill(randx[index], val);
00671    }
00672 
00673    // normalize
00674    hist->Scale( 1./hist->Integral("width") );
00675 
00676    return hist;
00677 }
00678 
00679 // *********************************************
00680 TH2D * BCIntegrate::MarginalizeByIntegrate(BCParameter * parameter1, BCParameter * parameter2)
00681 {
00682    // set parameter to marginalize
00683    ResetVarlist(1);
00684    int index1 = parameter1->GetIndex();
00685    UnsetVar(index1);
00686    int index2 = parameter2->GetIndex();
00687    UnsetVar(index2);
00688 
00689    // define histogram
00690    double hmin1 = parameter1->GetLowerLimit();
00691    double hmax1 = parameter1->GetUpperLimit();
00692    double hdx1  = (hmax1 - hmin1) / double(fNbins);
00693 
00694    double hmin2 = parameter2->GetLowerLimit();
00695    double hmax2 = parameter2->GetUpperLimit();
00696    double hdx2  = (hmax2 - hmin2) / double(fNbins);
00697 
00698    TH2D * hist = new TH2D(Form("hist_%s_%s", parameter1->GetName().data(), parameter2->GetName().data()),"",
00699          fNbins, hmin1, hmax1,
00700          fNbins, hmin2, hmax2);
00701 
00702    // integrate
00703    std::vector <double> randx;
00704    randx.assign(fNvar, 0.0);
00705 
00706    for(int i=0;i<=fNbins;i++) {
00707       randx[index1] = hmin1 + (double)i * hdx1;
00708       for(int j=0;j<=fNbins;j++) {
00709          randx[index2] = hmin2 + (double)j * hdx2;
00710 
00711          double val = IntegralMC(randx);
00712          hist->Fill(randx[index1],randx[index2], val);
00713       }
00714    }
00715 
00716    // normalize
00717    hist->Scale(1.0/hist->Integral("width"));
00718 
00719    return hist;
00720 }
00721 
00722 // *********************************************
00723 TH1D * BCIntegrate::MarginalizeByMetro(BCParameter * parameter)
00724 {
00725    int niter = fMarkovChainNIterations;
00726 
00727    if (fMarkovChainAutoN == true)
00728       niter = fNbins*fNbins*fNSamplesPer2DBin*fNvar;
00729 
00730    BCLog::OutDetail(Form(" --> Number of samples in Metropolis marginalization: %d.", niter));
00731 
00732    // set parameter to marginalize
00733    int index = parameter->GetIndex();
00734 
00735    // define histogram
00736    double hmin = parameter->GetLowerLimit();
00737    double hmax = parameter->GetUpperLimit();
00738    TH1D * hist = new TH1D("hist","", fNbins, hmin, hmax);
00739 
00740    // prepare Metro
00741    std::vector <double> randx;
00742    randx.assign(fNvar, 0.0);
00743    InitMetro();
00744 
00745    for(int i=0;i<=niter;i++) {
00746       GetRandomPointMetro(randx);
00747       hist->Fill(randx[index]);
00748    }
00749 
00750    // normalize
00751    hist->Scale(1.0/hist->Integral("width"));
00752 
00753    return hist;
00754 }
00755 
00756 // *********************************************
00757 TH2D * BCIntegrate::MarginalizeByMetro(BCParameter * parameter1, BCParameter * parameter2)
00758 {
00759    int niter=fNbins*fNbins*fNSamplesPer2DBin*fNvar;
00760 
00761    // set parameter to marginalize
00762    int index1 = parameter1->GetIndex();
00763    int index2 = parameter2->GetIndex();
00764 
00765    // define histogram
00766    double hmin1 = parameter1->GetLowerLimit();
00767    double hmax1 = parameter1->GetUpperLimit();
00768 
00769    double hmin2 = parameter2->GetLowerLimit();
00770    double hmax2 = parameter2->GetUpperLimit();
00771 
00772    TH2D * hist = new TH2D(Form("hist_%s_%s", parameter1->GetName().data(), parameter2->GetName().data()),"",
00773          fNbins, hmin1, hmax1,
00774          fNbins, hmin2, hmax2);
00775 
00776    // prepare Metro
00777    std::vector <double> randx;
00778    randx.assign(fNvar, 0.0);
00779    InitMetro();
00780 
00781    for(int i=0;i<=niter;i++) {
00782       GetRandomPointMetro(randx);
00783       hist->Fill(randx[index1],randx[index2]);
00784    }
00785 
00786    // normalize
00787    hist->Scale(1.0/hist->Integral("width"));
00788 
00789    return hist;
00790 }
00791 
00792 // *********************************************
00793 int BCIntegrate::MarginalizeAllByMetro(const char * name="")
00794 {
00795    int niter=fNbins*fNbins*fNSamplesPer2DBin*fNvar;
00796 
00797    BCLog::OutDetail(Form(" --> Number of samples in Metropolis marginalization: %d.", niter));
00798 
00799    // define 1D histograms
00800    for(int i=0;i<fNvar;i++)
00801    {
00802       double hmin1 = fx->at(i)->GetLowerLimit();
00803       double hmax1 = fx->at(i)->GetUpperLimit();
00804 
00805       TH1D * h1 = new TH1D(
00806             TString::Format("h%s_%s", name, fx->at(i)->GetName().data()),"",
00807             fNbins, hmin1, hmax1);
00808 
00809       fHProb1D.push_back(h1);
00810    }
00811 
00812    // define 2D histograms
00813    for(int i=0;i<fNvar-1;i++)
00814       for(int j=i+1;j<fNvar;j++) {
00815          double hmin1 = fx->at(i)->GetLowerLimit();
00816          double hmax1 = fx->at(i)->GetUpperLimit();
00817 
00818          double hmin2 = fx->at(j)->GetLowerLimit();
00819          double hmax2 = fx->at(j)->GetUpperLimit();
00820 
00821          TH2D * h2 = new TH2D(
00822             TString::Format("h%s_%s_%s", name, fx->at(i)->GetName().data(), fx->at(j)->GetName().data()),"",
00823             fNbins, hmin1, hmax1,
00824             fNbins, hmin2, hmax2);
00825 
00826          fHProb2D.push_back(h2);
00827       }
00828 
00829    // get number of 2d distributions
00830    int nh2d = fHProb2D.size();
00831 
00832    BCLog::OutDetail(Form(" --> Marginalizing %d 1D distributions and %d 2D distributions.", fNvar, nh2d));
00833 
00834    // prepare function fitting
00835    double dx = 0.;
00836    double dy = 0.;
00837 
00838    if (fFitFunctionIndexX >= 0) {
00839       dx = (fDataPointUpperBoundaries->GetValue(fFitFunctionIndexX) -
00840             fDataPointLowerBoundaries->GetValue(fFitFunctionIndexX))
00841             / double(fErrorBandNbinsX);
00842 
00843       dx = (fDataPointUpperBoundaries->GetValue(fFitFunctionIndexY) -
00844             fDataPointLowerBoundaries->GetValue(fFitFunctionIndexY))
00845             / double(fErrorBandNbinsY);
00846 
00847       fErrorBandXY = new TH2D(
00848             TString::Format("errorbandxy_%d",BCLog::GetHIndex()), "",
00849             fErrorBandNbinsX,
00850             fDataPointLowerBoundaries->GetValue(fFitFunctionIndexX) - 0.5 * dx,
00851             fDataPointUpperBoundaries->GetValue(fFitFunctionIndexX) + 0.5 * dx,
00852             fErrorBandNbinsY,
00853             fDataPointLowerBoundaries->GetValue(fFitFunctionIndexY) - 0.5 * dy,
00854             fDataPointUpperBoundaries->GetValue(fFitFunctionIndexY) + 0.5 * dy);
00855       fErrorBandXY->SetStats(kFALSE);
00856 
00857       for (int ix = 1; ix <= fErrorBandNbinsX; ++ix)
00858          for (int iy = 1; iy <= fErrorBandNbinsX; ++iy)
00859             fErrorBandXY->SetBinContent(ix, iy, 0.0);
00860    }
00861 
00862    // prepare Metro
00863    std::vector <double> randx;
00864 
00865    randx.assign(fNvar, 0.0);
00866    InitMetro();
00867 
00868    // reset counter
00869    fNacceptedMCMC = 0;
00870 
00871    // run Metro and fill histograms
00872    for(int i=0;i<=niter;i++) {
00873       GetRandomPointMetro(randx);
00874 
00875       // save this point to the markov chain in the ROOT file
00876       if (fFlagWriteMarkovChain) {
00877          fMCMCIteration = i;
00878          fMarkovChainTree->Fill();
00879       }
00880 
00881       for(int j=0;j<fNvar;j++)
00882          fHProb1D[j]->Fill( randx[j] );
00883 
00884       int ih=0;
00885       for(int j=0;j<fNvar-1;j++)
00886          for(int k=j+1;k<fNvar;k++) {
00887             fHProb2D[ih]->Fill(randx[j],randx[k]);
00888             ih++;
00889          }
00890 
00891       if((i+1)%100000==0)
00892          BCLog::OutDebug(Form("BCIntegrate::MarginalizeAllByMetro. %d samples done.",i+1));
00893 
00894       // function fitting
00895 
00896       if (fFitFunctionIndexX >= 0) {
00897          // loop over all possible x values ...
00898          if (fErrorBandContinuous) {
00899             double x = 0;
00900 
00901             for (int ix = 0; ix < fErrorBandNbinsX; ix++) {
00902                // calculate x
00903                x = fErrorBandXY->GetXaxis()->GetBinCenter(ix + 1);
00904 
00905                // calculate y
00906                std::vector <double> xvec;
00907                xvec.push_back(x);
00908                double y = FitFunction(xvec, randx);
00909                xvec.clear();
00910 
00911                // fill histogram
00912                fErrorBandXY->Fill(x, y);
00913             }
00914          }
00915 
00916          // ... or evaluate at the data point x-values
00917          else {
00918             int ndatapoints = int(fErrorBandX.size());
00919             double x = 0;
00920 
00921             for (int ix = 0; ix < ndatapoints; ++ix) {
00922                // calculate x
00923                x = fErrorBandX.at(ix);
00924 
00925                // calculate y
00926                std::vector <double> xvec;
00927                xvec.push_back(x);
00928                double y = FitFunction(xvec, randx);
00929                xvec.clear();
00930 
00931                // fill histogram
00932                fErrorBandXY->Fill(x, y);
00933             }
00934          }
00935       }
00936    }
00937 
00938    // normalize histograms
00939    for(int i=0;i<fNvar;i++)
00940       fHProb1D[i]->Scale( 1./fHProb1D[i]->Integral("width") );
00941 
00942    for (int i=0;i<nh2d;i++)
00943       fHProb2D[i]->Scale( 1./fHProb2D[i]->Integral("width") );
00944 
00945    if (fFitFunctionIndexX >= 0)
00946       fErrorBandXY->Scale(1.0/fErrorBandXY->Integral() * fErrorBandXY->GetNbinsX());
00947 
00948    if (fFitFunctionIndexX >= 0) {
00949       for (int ix = 1; ix <= fErrorBandNbinsX; ix++) {
00950          double sum = 0;
00951 
00952          for (int iy = 1; iy <= fErrorBandNbinsY; iy++)
00953             sum += fErrorBandXY->GetBinContent(ix, iy);
00954 
00955          for (int iy = 1; iy <= fErrorBandNbinsY; iy++) {
00956             double newvalue = fErrorBandXY->GetBinContent(ix, iy) / sum;
00957             fErrorBandXY->SetBinContent(ix, iy, newvalue);
00958          }
00959       }
00960    }
00961 
00962    BCLog::OutDebug(Form("BCIntegrate::MarginalizeAllByMetro done with %i trials and %i accepted trials. Efficiency is %f",fNmetro, fNacceptedMCMC, double(fNacceptedMCMC)/double(fNmetro)));
00963 
00964    return fNvar+nh2d;
00965 }
00966 
00967 // *********************************************
00968 TH1D * BCIntegrate::GetH1D(int parIndex)
00969 {
00970    if(fHProb1D.size()==0) {
00971       BCLog::OutWarning("BCModel::GetH1D. MarginalizeAll() has to be run prior to this to fill the distributions.");
00972       return 0;
00973    }
00974 
00975    if(parIndex<0 || parIndex>fNvar) {
00976       BCLog::OutWarning(Form("BCIntegrate::GetH1D. Parameter index %d is invalid.",parIndex));
00977       return 0;
00978    }
00979 
00980    return fHProb1D[parIndex];
00981 }
00982 
00983 // *********************************************
00984 int BCIntegrate::GetH2DIndex(int parIndex1, int parIndex2)
00985 {
00986    if(parIndex1>fNvar-1 || parIndex1<0) {
00987       BCLog::OutWarning(Form("BCIntegrate::GetH2DIndex. Parameter index %d is invalid", parIndex1));
00988       return -1;
00989    }
00990 
00991    if(parIndex2>fNvar-1 || parIndex2<0) {
00992       BCLog::OutWarning(Form("BCIntegrate::GetH2DIndex. Parameter index %d is invalid", parIndex2));
00993       return -1;
00994    }
00995 
00996    if(parIndex1==parIndex2) {
00997       BCLog::OutWarning(Form("BCIntegrate::GetH2DIndex. Parameters have equal indeces: %d , %d", parIndex1, parIndex2));
00998       return -1;
00999    }
01000 
01001    if(parIndex1>parIndex2) {
01002       BCLog::OutWarning("BCIntegrate::GetH2DIndex. First parameters must be smaller than second (sorry).");
01003       return -1;
01004    }
01005 
01006    int index=0;
01007    for(int i=0;i<fNvar-1;i++)
01008       for(int j=i+1;j<fNvar;j++) {
01009          if(i==parIndex1 && j==parIndex2)
01010             return index;
01011          index++;
01012       }
01013 
01014    BCLog::OutWarning("BCIntegrate::GetH2DIndex. Invalid index combination.");
01015 
01016    return -1;
01017 }
01018 
01019 // *********************************************
01020 TH2D * BCIntegrate::GetH2D(int parIndex1, int parIndex2)
01021 {
01022    if(fHProb2D.size()==0) {
01023       BCLog::OutWarning("BCModel::GetH2D. MarginalizeAll() has to be run prior to this to fill the distributions.");
01024       return 0;
01025    }
01026 
01027    int hindex = GetH2DIndex(parIndex1, parIndex2);
01028    if(hindex==-1)
01029       return 0;
01030 
01031    if(hindex>(int)fHProb2D.size()-1) {
01032       BCLog::OutWarning("BCIntegrate::GetH2D. Got invalid index from GetH2DIndex(). Something went wrong.");
01033       return 0;
01034    }
01035 
01036    return fHProb2D[hindex];
01037 }
01038 
01039 // *********************************************
01040 double BCIntegrate::GetRandomPoint(std::vector <double> &x)
01041 {
01042    GetRandomVector(x);
01043 
01044    for(int i=0;i<fNvar;i++)
01045       x[i]=fMin[i]+x[i]*(fMax[i]-fMin[i]);
01046 
01047    return Eval(x);
01048 }
01049 
01050 // *********************************************
01051 double BCIntegrate::GetRandomPointImportance(std::vector <double> &x)
01052 {
01053    double p = 1.1;
01054    double g = 1.0;
01055 
01056    while (p>g) {
01057       GetRandomVector(x);
01058 
01059       for(int i=0;i<fNvar;i++)
01060          x[i] = fMin[i] + x[i] * (fMax[i]-fMin[i]);
01061 
01062       p = fRandom->Rndm();
01063 
01064       g = EvalSampling(x);
01065    }
01066 
01067    return Eval(x);
01068 }
01069 
01070 // *********************************************
01071 void BCIntegrate::InitMetro()
01072 {
01073    fNmetro=0;
01074 
01075    if (fXmetro0.size() <= 0) {
01076       // start in the center of the phase space
01077       for(int i=0;i<fNvar;i++)
01078          fXmetro0.push_back((fMin[i]+fMax[i])/2.0);
01079    }
01080 
01081    fMarkovChainValue =  LogEval(fXmetro0);
01082 
01083    // run metropolis for a few times and dump the result... (to forget the initial position)
01084    std::vector <double> x;
01085    x.assign(fNvar,0.);
01086 
01087    for(int i=0;i<1000;i++)
01088       GetRandomPointMetro(x);
01089 
01090    fNmetro = 0;
01091 }
01092 
01093 // *********************************************
01094 void BCIntegrate::GetRandomPointMetro(std::vector <double> &x)
01095 {
01096    // get new point
01097    GetRandomVectorMetro(fXmetro1);
01098 
01099    // scale the point to the allowed region and stepsize
01100    int in=1;
01101    for(int i=0;i<fNvar;i++) {
01102       fXmetro1[i] = fXmetro0[i] + fXmetro1[i] * (fMax[i]-fMin[i]);
01103 
01104       // check whether the generated point is inside the allowed region
01105       if( fXmetro1[i]<fMin[i] || fXmetro1[i]>fMax[i] )
01106          in=0;
01107    }
01108 
01109    // calculate the log probabilities and compare old and new point
01110 
01111    double p0 = fMarkovChainValue; // old point
01112    double p1 = 0; // new point
01113    int accept=0;
01114 
01115    if(in) {
01116       p1 = LogEval(fXmetro1);
01117 
01118       if(p1>=p0)
01119          accept=1;
01120       else {
01121          double r=log(fRandom->Rndm());
01122          if(r<p1-p0)
01123             accept=1;
01124       }
01125    }
01126 
01127    // fill the return point after the decision
01128    if(accept) {
01129       // increase counter
01130       fNacceptedMCMC++;
01131 
01132       for(int i=0;i<fNvar;i++) {
01133          fXmetro0[i]=fXmetro1[i];
01134          x[i]=fXmetro1[i];
01135          fMarkovChainValue = p1;
01136       }
01137    }
01138    else
01139       for(int i=0;i<fNvar;i++) {
01140          x[i]=fXmetro0[i];
01141          fXmetro1[i] = x[i];
01142          fMarkovChainValue = p0;
01143       }
01144 
01145    fNmetro++;
01146 }
01147 
01148 // *********************************************
01149 void BCIntegrate::GetRandomPointSamplingMetro(std::vector <double> &x)
01150 {
01151    // get new point
01152    GetRandomVectorMetro(fXmetro1);
01153 
01154    // scale the point to the allowed region and stepsize
01155    int in=1;
01156    for(int i=0;i<fNvar;i++) {
01157       fXmetro1[i] = fXmetro0[i] + fXmetro1[i] * (fMax[i]-fMin[i]);
01158 
01159       // check whether the generated point is inside the allowed region
01160       if( fXmetro1[i]<fMin[i] || fXmetro1[i]>fMax[i] )
01161          in=0;
01162    }
01163 
01164    // calculate the log probabilities and compare old and new point
01165    double p0 = LogEvalSampling(fXmetro0); // old point
01166    double p1=0; // new point
01167    if(in)
01168       p1= LogEvalSampling(fXmetro1);
01169 
01170    // compare log probabilities
01171    int accept=0;
01172    if(in) {
01173       if(p1>=p0)
01174          accept=1;
01175       else {
01176          double r=log(fRandom->Rndm());
01177          if(r<p1-p0)
01178             accept=1;
01179       }
01180    }
01181 
01182    // fill the return point after the decision
01183    if(accept)
01184       for(int i=0;i<fNvar;i++) {
01185          fXmetro0[i]=fXmetro1[i];
01186          x[i]=fXmetro1[i];
01187       }
01188    else
01189       for(int i=0;i<fNvar;i++)
01190          x[i]=fXmetro0[i];
01191 
01192    fNmetro++;
01193 }
01194 
01195 // *********************************************
01196 void BCIntegrate::GetRandomVector(std::vector <double> &x)
01197 {
01198    double * randx = new double[fNvar];
01199 
01200    fRandom->RndmArray(fNvar, randx);
01201 
01202    for(int i=0;i<fNvar;i++)
01203       x[i] = randx[i];
01204 
01205    delete[] randx;
01206    randx = 0;
01207 }
01208 
01209 // *********************************************
01210 void BCIntegrate::GetRandomVectorMetro(std::vector <double> &x)
01211 {
01212    double * randx = new double[fNvar];
01213 
01214    fRandom->RndmArray(fNvar, randx);
01215 
01216    for(int i=0;i<fNvar;i++)
01217       x[i] = 2.0 * (randx[i] - 0.5) * fMarkovChainStepSize;
01218 
01219    delete[] randx;
01220    randx = 0;
01221 }
01222 
01223 // *********************************************
01224 TMinuit * BCIntegrate::GetMinuit()
01225 {
01226    if (!fMinuit)
01227       fMinuit = new TMinuit();
01228 
01229    return fMinuit;
01230 }
01231 
01232 // *********************************************
01233 
01234 void BCIntegrate::FindModeMinuit(std::vector<double> start, int printlevel)
01235 {
01236    bool have_start = true;
01237 
01238    // check start values
01239    if (int(start.size()) != fNvar)
01240       have_start = false;
01241 
01242    // set global this
01243    global_this = this;
01244 
01245    // define minuit
01246    if (fMinuit)
01247       delete fMinuit;
01248    fMinuit = new TMinuit(fNvar);
01249 
01250    // set print level
01251    fMinuit->SetPrintLevel(printlevel);
01252 
01253    // set function
01254    fMinuit->SetFCN(&BCIntegrate::FCNLikelihood);
01255 
01256    // set UP for likelihood
01257    fMinuit->SetErrorDef(0.5);
01258 
01259    // set parameters
01260    int flag;
01261    for (int i = 0; i < fNvar; i++) {
01262       double starting_point = (fMin[i]+fMax[i])/2.;
01263       if(have_start)
01264          starting_point = start[i];
01265          fMinuit->mnparm(i,
01266                          fx->at(i)->GetName().data(),
01267                          starting_point,
01268                          (fMax[i]-fMin[i])/100.,
01269                          fMin[i],
01270                          fMax[i],
01271                          flag);
01272    }
01273 
01274    // do mcmc minimization
01275    //   fMinuit->mnseek();
01276 
01277    // do minimization
01278    fMinuit->mnexcm("MIGRAD", fMinuitArglist, 2, flag);
01279 
01280    // improve search for local minimum
01281    //   fMinuit->mnimpr();
01282 
01283    // copy flag
01284    fMinuitErrorFlag = flag;
01285 
01286    // check if mode found by minuit is better than previous estimation
01287    double probmax = 0;
01288    bool valid = false;
01289 
01290    if ( int(fBestFitParameters.size()) == fNvar) {
01291       valid = true;
01292       for (int i = 0; i < fNvar; ++i)
01293          if (fBestFitParameters.at(i) < fMin[i] || fBestFitParameters.at(i) > fMax[i])
01294             valid= false;
01295 
01296       if (valid)
01297          probmax = Eval( fBestFitParameters );
01298    }
01299 
01300    std::vector<double> tempvec;
01301    for (int i = 0; i < fNvar; i++) {
01302       double par;
01303       double parerr;
01304       fMinuit->GetParameter(i, par, parerr);
01305       tempvec.push_back(par);
01306    }
01307    double probmaxminuit = Eval( tempvec );
01308 
01309    // set best fit parameters
01310    if ( (probmaxminuit > probmax && !fFlagIgnorePrevOptimization) || !valid || fFlagIgnorePrevOptimization) {
01311       fBestFitParameters.clear();
01312       fBestFitParameterErrors.clear();
01313 
01314       for (int i = 0; i < fNvar; i++) {
01315          double par;
01316          double parerr;
01317          fMinuit->GetParameter(i, par, parerr);
01318          fBestFitParameters.push_back(par);
01319          fBestFitParameterErrors.push_back(parerr);
01320       }
01321 
01322       // set optimization method used to find the mode
01323       fOptimizationMethodMode = BCIntegrate::kOptMinuit;
01324    }
01325 
01326    // delete minuit
01327    delete fMinuit;
01328    fMinuit = 0;
01329 
01330    return;
01331 }
01332 
01333 // *********************************************
01334 void BCIntegrate::InitializeSATree()
01335 {
01336    delete fTreeSA;
01337    fTreeSA = new TTree("SATree", "SATree");
01338 
01339    fTreeSA->Branch("Iteration",      &fSANIterations,   "iteration/I");
01340    fTreeSA->Branch("NParameters",    &fNvar,            "parameters/I");
01341    fTreeSA->Branch("Temperature",    &fSATemperature,   "temperature/D");
01342    fTreeSA->Branch("LogProbability", &fSALogProb,       "log(probability)/D");
01343 
01344    for (int i = 0; i < fNvar; ++i)
01345       fTreeSA->Branch(TString::Format("Parameter%i", i), &fSAx[i], TString::Format("parameter %i/D", i));
01346 }
01347 
01348 // *********************************************
01349 void BCIntegrate::FindModeSA(std::vector<double> start)
01350 {
01351    // note: if f(x) is the function to be minimized, then
01352    // f(x) := - LogEval(parameters)
01353 
01354    bool have_start = true;
01355    std::vector<double> x, y, best_fit; // vectors for current state, new proposed state and best fit up to now
01356    double fval_x, fval_y, fval_best_fit; // function values at points x, y and best_fit (we save them rather than to re-calculate them every time)
01357    int t = 1; // time iterator
01358 
01359    // check start values
01360    if (int(start.size()) != fNvar)
01361       have_start = false;
01362 
01363    // if no starting point is given, set to center of parameter space
01364    if ( !have_start ) {
01365       start.clear();
01366       for (int i = 0; i < fNvar; i++)
01367          start.push_back((fMin[i]+fMax[i])/2.);
01368    }
01369 
01370    // set current state and best fit to starting point
01371    x.clear();
01372    best_fit.clear();
01373    for (int i = 0; i < fNvar; i++) {
01374       x.push_back(start[i]);
01375       best_fit.push_back(start[i]);
01376    }
01377    // calculate function value at starting point
01378    fval_x = fval_best_fit = LogEval(x);
01379 
01380    // run while still "hot enough"
01381    while ( SATemperature(t) > fSATmin ) {
01382       // generate new state
01383       y = GetProposalPointSA(x, t);
01384 
01385       // check if the proposed point is inside the phase space
01386       // if not, reject it
01387       bool is_in_ranges = true;
01388 
01389       for (int i = 0; i < fNvar; i++)
01390          if (y[i] > fMax[i] || y[i] < fMin[i])
01391             is_in_ranges = false;
01392 
01393       if ( !is_in_ranges )
01394          ; // do nothing...
01395       else {
01396          // calculate function value at new point
01397          fval_y = LogEval(y);
01398 
01399          // is it better than the last one?
01400          // if so, update state and chef if it is the new best fit...
01401          if (fval_y >= fval_x) {
01402             x.clear();
01403             for (int i = 0; i < fNvar; i++)
01404                x.push_back(y[i]);
01405 
01406             fval_x = fval_y;
01407 
01408             if (fval_y > fval_best_fit) {
01409                best_fit.clear();
01410                for (int i = 0; i < fNvar; i++)
01411                   best_fit.push_back(y[i]);
01412 
01413                fval_best_fit = fval_y;
01414             }
01415          }
01416          // ...else, only accept new state w/ certain probability
01417          else {
01418             if (fRandom->Rndm() <= exp( (fval_y - fval_x) / SATemperature(t) )) {
01419                x.clear();
01420                for (int i = 0; i < fNvar; i++)
01421                   x.push_back(y[i]);
01422 
01423                fval_x = fval_y;
01424             }
01425          }
01426       }
01427 
01428       // update tree variables
01429       fSANIterations = t;
01430       fSATemperature = SATemperature(t);
01431       fSALogProb = fval_x;
01432       fSAx.clear();
01433       for (int i = 0; i < fNvar; ++i)
01434          fSAx.push_back(x[i]);
01435 
01436       // fill tree
01437       if (fFlagWriteSAToFile)
01438          fTreeSA->Fill();
01439 
01440       // increate t
01441       t++;
01442    }
01443 
01444    // check if mode found by minuit is better than previous estimation
01445    double probmax = 0;
01446    bool valid = false;
01447 
01448    if ( int(fBestFitParameters.size()) == fNvar) {
01449       valid = true;
01450       for (int i = 0; i < fNvar; ++i)
01451          if (fBestFitParameters.at(i) < fMin[i] || fBestFitParameters.at(i) > fMax[i])
01452             valid= false;
01453 
01454       if (valid)
01455          probmax = Eval( fBestFitParameters );
01456    }
01457 
01458    double probmaxsa = Eval( best_fit);
01459 
01460    if ( (probmaxsa > probmax  && !fFlagIgnorePrevOptimization) || !valid || fFlagIgnorePrevOptimization) {
01461       // set best fit parameters
01462       fBestFitParameters.clear();
01463       fBestFitParameterErrors.clear();
01464 
01465       for (int i = 0; i < fNvar; i++) {
01466          fBestFitParameters.push_back(best_fit[i]);
01467          fBestFitParameterErrors.push_back(-1.); // error undefined
01468       }
01469 
01470       // set optimization moethod used to find the mode
01471       fOptimizationMethodMode = BCIntegrate::kOptSA;
01472    }
01473 
01474    return;
01475 }
01476 
01477 // *********************************************
01478 
01479 double BCIntegrate::SATemperature(double t)
01480 {
01481    // do we have Cauchy (default), Boltzmann or custom annealing schedule?
01482    if (fSASchedule == BCIntegrate::kSABoltzmann)
01483       return SATemperatureBoltzmann(t);
01484    else if (fSASchedule == BCIntegrate::kSACauchy)
01485       return SATemperatureCauchy(t);
01486    else
01487       return SATemperatureCustom(t);
01488 }
01489 
01490 // *********************************************
01491 
01492 double BCIntegrate::SATemperatureBoltzmann(double t)
01493 {
01494    return fSAT0 / log((double)(t + 1));
01495 }
01496 
01497 // *********************************************
01498 
01499 double BCIntegrate::SATemperatureCauchy(double t)
01500 {
01501    return fSAT0 / (double)t;
01502 }
01503 
01504 // *********************************************
01505 
01506 double BCIntegrate::SATemperatureCustom(double t)
01507 {
01508    BCLog::OutError("BCIntegrate::SATemperatureCustom : No custom temperature schedule defined");
01509    return 0.;
01510 }
01511 
01512 // *********************************************
01513 
01514 std::vector<double> BCIntegrate::GetProposalPointSA(std::vector<double> x, int t)
01515 {
01516    // do we have Cauchy (default), Boltzmann or custom annealing schedule?
01517    if (fSASchedule == BCIntegrate::kSABoltzmann)
01518       return GetProposalPointSABoltzmann(x, t);
01519    else if (fSASchedule == BCIntegrate::kSACauchy)
01520       return GetProposalPointSACauchy(x, t);
01521    else
01522       return GetProposalPointSACustom(x, t);
01523 }
01524 
01525 // *********************************************
01526 
01527 std::vector<double> BCIntegrate::GetProposalPointSABoltzmann(std::vector<double> x, int t)
01528 {
01529    std::vector<double> y;
01530    y.clear();
01531    double new_val, norm;
01532 
01533    for (int i = 0; i < fNvar; i++) {
01534       norm = (fMax[i] - fMin[i]) * SATemperature(t) / 2.;
01535       new_val = x[i] + norm * fRandom->Gaus();
01536       y.push_back(new_val);
01537    }
01538 
01539    return y;
01540 }
01541 
01542 // *********************************************
01543 
01544 std::vector<double> BCIntegrate::GetProposalPointSACauchy(std::vector<double> x, int t)
01545 {
01546    std::vector<double> y;
01547    y.clear();
01548 
01549    if (fNvar == 1) {
01550       double cauchy, new_val, norm;
01551 
01552       norm = (fMax[0] - fMin[0]) * SATemperature(t) / 2.;
01553       cauchy = tan(3.14159 * (fRandom->Rndm() - 0.5));
01554       new_val = x[0] + norm * cauchy;
01555       y.push_back(new_val);
01556    }
01557    else {
01558       // use sampling to get radial n-dim Cauchy distribution
01559 
01560       // first generate a random point uniformly distributed on a
01561       // fNvar-dimensional hypersphere
01562       y = SAHelperGetRandomPointOnHypersphere();
01563 
01564       // scale the vector by a random factor determined by the radial
01565       // part of the fNvar-dimensional Cauchy distribution
01566       double radial = SATemperature(t) * SAHelperGetRadialCauchy();
01567 
01568       // scale y by radial part and the size of dimension i in phase space
01569       // afterwards, move by x
01570       for (int i = 0; i < fNvar; i++)
01571          y[i] = (fMax[i] - fMin[i]) * y[i] * radial / 2. + x[i];
01572    }
01573 
01574    return y;
01575 }
01576 
01577 // *********************************************
01578 
01579 std::vector<double> BCIntegrate::GetProposalPointSACustom(std::vector<double> x, int t)
01580 {
01581    BCLog::OutError("BCIntegrate::GetProposalPointSACustom : No custom proposal function defined");
01582    return std::vector<double>(fNvar);
01583 }
01584 
01585 // *********************************************
01586 
01587 std::vector<double> BCIntegrate::SAHelperGetRandomPointOnHypersphere()
01588 {
01589    std::vector<double> rand_point(fNvar);
01590 
01591    // This method can only be called with fNvar >= 2 since the 1-dim case
01592    // is already hard wired into the Cauchy annealing proposal function.
01593    // To speed things up, hard-code fast method for 2 and dimensions.
01594    // The algorithm for 2D can be found at
01595    // http://mathworld.wolfram.com/CirclePointPicking.html
01596    // For 3D just using ROOT's algorithm.
01597    if (fNvar == 2) {
01598       double x1, x2, s;
01599       do {
01600          x1 = fRandom->Rndm() * 2. - 1.;
01601          x2 = fRandom->Rndm() * 2. - 1.;
01602          s = x1*x1 + x2*x2;
01603       }
01604       while (s >= 1);
01605 
01606       rand_point[0] = (x1*x1 - x2*x2) / s;
01607       rand_point[1] = (2.*x1*x2) / s;
01608    }
01609    else if (fNvar == 3) {
01610       fRandom->Sphere(rand_point[0], rand_point[1], rand_point[2], 1.0);
01611    }
01612    else {
01613       double s = 0.,
01614          gauss_num;
01615 
01616       for (int i = 0; i < fNvar; i++) {
01617          gauss_num = fRandom->Gaus();
01618          rand_point[i] = gauss_num;
01619          s += gauss_num * gauss_num;
01620       }
01621       s = sqrt(s);
01622 
01623       for (int i = 0; i < fNvar; i++)
01624          rand_point[i] = rand_point[i] / s;
01625    }
01626 
01627    return rand_point;
01628 }
01629 
01630 // *********************************************
01631 
01632 double BCIntegrate::SAHelperGetRadialCauchy()
01633 {
01634    // theta is sampled from a rather complicated distribution,
01635    // so first we create a lookup table with 10000 random numbers
01636    // once and then, each time we need a new random number,
01637    // we just look it up in the table.
01638    double theta;
01639 
01640    // static vectors for theta-sampling-map
01641    static double map_u[10001];
01642    static double map_theta[10001];
01643    static bool initialized = false;
01644    static int map_dimension = 0;
01645 
01646    // is the lookup-table already initialized? if not, do it!
01647    if (!initialized || map_dimension != fNvar) {
01648       double init_theta;
01649       double beta = SAHelperSinusToNIntegral(fNvar - 1, 1.57079632679);
01650 
01651       for (int i = 0; i <= 10000; i++) {
01652          init_theta = 3.14159265 * (double)i / 5000.;
01653          map_theta[i] = init_theta;
01654 
01655          map_u[i] = SAHelperSinusToNIntegral(fNvar - 1, init_theta) / beta;
01656       }
01657 
01658       map_dimension = fNvar;
01659       initialized = true;
01660    } // initializing is done.
01661 
01662    // generate uniform random number for sampling
01663    double u = fRandom->Rndm();
01664 
01665    // Find the two elements just greater than and less than u
01666    // using a binary search (O(log(N))).
01667    int lo = 0;
01668    int up = 10000;
01669    int mid;
01670 
01671    while (up != lo) {
01672       mid = ((up - lo + 1) / 2) + lo;
01673 
01674       if (u >= map_u[mid])
01675          lo = mid;
01676       else
01677          up = mid - 1;
01678    }
01679    up++;
01680 
01681    // perform linear interpolation:
01682    theta = map_theta[lo] + (u - map_u[lo]) / (map_u[up] - map_u[lo]) * (map_theta[up] - map_theta[lo]);
01683 
01684    return tan(theta);
01685 }
01686 
01687 // *********************************************
01688 
01689 double BCIntegrate::SAHelperSinusToNIntegral(int dim, double theta)
01690 {
01691    if (dim < 1)
01692       return theta;
01693    else if (dim == 1)
01694       return (1. - cos(theta));
01695    else if (dim == 2)
01696       return 0.5 * (theta - sin(theta) * cos(theta));
01697    else if (dim == 3)
01698       return (2. - sin(theta) * sin(theta) * cos(theta) - 2. * cos(theta)) / 3.;
01699    else
01700       return - pow(sin(theta), (double)(dim - 1)) * cos(theta) / (double)dim
01701          + (double)(dim - 1) / (double)dim
01702          * SAHelperSinusToNIntegral(dim - 2, theta);
01703 }
01704 
01705 // *********************************************
01706 
01707 void BCIntegrate::SetMode(std::vector <double> mode)
01708 {
01709    if((int)mode.size() == fNvar) {
01710       fBestFitParameters.clear();
01711       for (int i = 0; i < fNvar; i++)
01712          fBestFitParameters.push_back(mode[i]);
01713    }
01714 }
01715 
01716 // *********************************************
01717 
01718 void BCIntegrate::FCNLikelihood(int &npar, double * grad, double &fval, double * par, int flag)
01719 {
01720    // copy parameters
01721    std::vector <double> parameters;
01722 
01723    int n = global_this->GetNvar();
01724 
01725    for (int i = 0; i < n; i++)
01726       parameters.push_back(par[i]);
01727 
01728    fval = - global_this->LogEval(parameters);
01729 
01730    // delete parameters
01731    parameters.clear();
01732 }
01733 
01734 // *********************************************
01735 
01736 //void BCIntegrate::FindModeMCMC(int flag_run)
01737 void BCIntegrate::FindModeMCMC()
01738 {
01739    // call PreRun
01740 //   if (flag_run == 0)
01741    if (!fMCMCFlagPreRun)
01742       MCMCMetropolisPreRun();
01743 
01744    // find global maximum
01745    //   double probmax = (MCMCGetMaximumLogProb()).at(0);
01746 
01747    double probmax = 0;
01748 
01749    if ( int(fBestFitParameters.size()) == fNvar) {
01750       probmax = Eval( fBestFitParameters );
01751 //      fBestFitParameters = MCMCGetMaximumPoint(0);
01752    }
01753 
01754    // loop over all chains and find the maximum point
01755    for (int i = 0; i < fMCMCNChains; ++i) {
01756       double prob = exp( (MCMCGetMaximumLogProb()).at(i));
01757 
01758       // copy the point into the vector
01759       if ( (prob >= probmax && !fFlagIgnorePrevOptimization) || fFlagIgnorePrevOptimization) {
01760          probmax = prob;
01761 
01762          fBestFitParameters.clear();
01763          fBestFitParameterErrors.clear();
01764          fBestFitParameters = MCMCGetMaximumPoint(i);
01765 
01766          for (int j = 0; j < fNvar; ++j)
01767             fBestFitParameterErrors.push_back(-1.); // error undefined
01768 
01769          fOptimizationMethodMode = BCIntegrate::kOptMetropolis;
01770       }
01771    }
01772 }
01773 
01774 // *********************************************
01775 void BCIntegrate::CubaIntegrand(const int *ndim, const double xx[],
01776       const int *ncomp, double ff[])
01777 {
01778 #ifdef HAVE_CUBA_H
01779    // scale variables
01780    double jacobian = 1.0;
01781 
01782    std::vector<double> scaled_parameters;
01783 
01784    for (int i = 0; i < *ndim; i++) {
01785       double range = global_this->fx->at(i)->GetUpperLimit() -  global_this->fx->at(i)->GetLowerLimit();
01786 
01787       // multiply range to jacobian
01788       jacobian *= range;
01789 
01790       // get the scaled parameter value
01791       scaled_parameters.push_back(global_this->fx->at(i)->GetLowerLimit() + xx[i] * range);
01792    }
01793 
01794    // call function to integrate
01795    ff[0] = global_this->Eval(scaled_parameters);
01796 
01797    // multiply jacobian
01798    ff[0] *= jacobian;
01799 
01800    // multiply fudge factor
01801    ff[0] *= 1e99;
01802 
01803    // remove parameter vector
01804    scaled_parameters.clear();
01805 #else
01806    BCLog::OutError("!!! This version of BAT is compiled without Cuba.");
01807    BCLog::OutError("    Use other integration methods or install Cuba and recompile BAT.");
01808 #endif
01809 }
01810 
01811 // *********************************************
01812 double BCIntegrate::CubaIntegrate()
01813 {
01814 #ifdef HAVE_CUBA_H
01815    double EPSREL = 1e-3;
01816    double EPSABS = 1e-12;
01817    double VERBOSE   = 0;
01818    double MINEVAL   = 0;
01819    double MAXEVAL   = 2000000;
01820    double NSTART    = 25000;
01821    double NINCREASE = 25000;
01822 
01823    std::vector<double> parameters_double;
01824    std::vector<double>    parameters_int;
01825 
01826    parameters_double.push_back(EPSREL);
01827    parameters_double.push_back(EPSABS);
01828 
01829    parameters_int.push_back(VERBOSE);
01830    parameters_int.push_back(MINEVAL);
01831    parameters_int.push_back(MAXEVAL);
01832    parameters_int.push_back(NSTART);
01833    parameters_int.push_back(NINCREASE);
01834 
01835    // print to log
01836    BCLog::OutDebug( Form(" --> Running Cuba/Vegas integation over %i dimensions.", fNvar));
01837    BCLog::OutDebug( Form(" --> Maximum number of iterations: %i", (int)MAXEVAL));
01838    BCLog::OutDebug( Form(" --> Aimed relative precision:     %e", EPSREL));
01839 
01840    return CubaIntegrate(0, parameters_double, parameters_int);
01841 #else
01842    BCLog::OutError("!!! This version of BAT is compiled without Cuba.");
01843    BCLog::OutError("    Use other integration methods or install Cuba and recompile BAT.");
01844    return 0.;
01845 #endif
01846 }
01847 
01848 // *********************************************
01849 double BCIntegrate::CubaIntegrate(int method, std::vector<double> parameters_double, std::vector<double> parameters_int)
01850 {
01851 #ifdef HAVE_CUBA_H
01852    const int NDIM      = int(fx ->size());
01853    const int NCOMP     = 1;
01854 
01855    const double EPSREL = parameters_double[0];
01856    const double EPSABS = parameters_double[1];
01857    const int VERBOSE   = int(parameters_int[0]);
01858    const int MINEVAL   = int(parameters_int[1]);
01859    const int MAXEVAL   = int(parameters_int[2]);
01860 
01861    int neval;
01862    int fail;
01863    int nregions;
01864    double integral[NCOMP];
01865    double error[NCOMP];
01866    double prob[NCOMP];
01867 
01868    global_this = this;
01869 
01870    integrand_t an_integrand = &BCIntegrate::CubaIntegrand;
01871 
01872    if (method == 0) {
01873       // set VEGAS specific parameters
01874       const int NSTART    = int(parameters_int[3]);
01875       const int NINCREASE = int(parameters_int[4]);
01876 
01877       // call VEGAS integration method
01878       Vegas(NDIM, NCOMP, an_integrand,
01879          EPSREL, EPSABS, VERBOSE, MINEVAL, MAXEVAL,
01880          NSTART, NINCREASE,
01881          &neval, &fail, integral, error, prob);
01882    }
01883    else if (method == 1) {
01884       // set SUAVE specific parameters
01885       const int LAST     = int(parameters_int[3]);
01886       const int NNEW     = int(parameters_int[4]);
01887       const int FLATNESS = int(parameters_int[5]);
01888 
01889       // call SUAVE integration method
01890       Suave(NDIM, NCOMP, an_integrand,
01891          EPSREL, EPSABS, VERBOSE | LAST, MINEVAL, MAXEVAL,
01892          NNEW, FLATNESS,
01893          &nregions, &neval, &fail, integral, error, prob);
01894    }
01895    else if (method == 2) {
01896       // set DIVONNE specific parameters
01897       const int KEY1         = int(parameters_int[3]);
01898       const int KEY2         = int(parameters_int[4]);
01899       const int KEY3         = int(parameters_int[5]);
01900       const int MAXPASS      = int(parameters_int[6]);
01901       const int BORDER       = int(parameters_int[7]);
01902       const int MAXCHISQ     = int(parameters_int[8]);
01903       const int MINDEVIATION = int(parameters_int[9]);
01904       const int NGIVEN       = int(parameters_int[10]);
01905       const int LDXGIVEN     = int(parameters_int[11]);
01906       const int NEXTRA       = int(parameters_int[12]);
01907 
01908       // call DIVONNE integration method
01909       Divonne(NDIM, NCOMP, an_integrand,
01910          EPSREL, EPSABS, VERBOSE, MINEVAL, MAXEVAL,
01911          KEY1, KEY2, KEY3, MAXPASS, BORDER, MAXCHISQ, MINDEVIATION,
01912          NGIVEN, LDXGIVEN, NULL, NEXTRA, NULL,
01913          &nregions, &neval, &fail, integral, error, prob);
01914    }
01915    else if (method == 3) {
01916       // set CUHRE specific parameters
01917       const int LAST = int(parameters_int[3]);
01918       const int KEY  = int(parameters_int[4]);
01919 
01920       // call CUHRE integration method
01921       Cuhre(NDIM, NCOMP, an_integrand,
01922          EPSREL, EPSABS, VERBOSE | LAST, MINEVAL, MAXEVAL, KEY,
01923          &nregions, &neval, &fail, integral, error, prob);
01924    }
01925    else {
01926       std::cout << " Integration method not available. " << std::endl;
01927       integral[0] = -1e99;
01928    }
01929 
01930    if (fail != 0) {
01931       std::cout << " Warning, integral did not converge with the given set of parameters. "<< std::endl;
01932       std::cout << " neval    = " << neval       << std::endl;
01933       std::cout << " fail     = " << fail        << std::endl;
01934       std::cout << " integral = " << integral[0] << std::endl;
01935       std::cout << " error    = " << error[0]    << std::endl;
01936       std::cout << " prob     = " << prob[0]     << std::endl;
01937    }
01938 
01939    return integral[0] / 1e99;
01940 #else
01941    BCLog::OutError("!!! This version of BAT is compiled without Cuba.");
01942    BCLog::OutError("    Use other integration methods or install Cuba and recompile BAT.");
01943    return 0.;
01944 #endif
01945 }
01946 
01947 // *********************************************
01948 void BCIntegrate::MCMCIterationInterface()
01949 {
01950    // what's within this method will be executed
01951    // for every iteration of the MCMC
01952 
01953    // fill error band
01954    MCMCFillErrorBand();
01955 
01956    // do user defined stuff
01957    MCMCUserIterationInterface();
01958 }
01959 
01960 // *********************************************
01961 void BCIntegrate::MCMCFillErrorBand()
01962 {
01963    if (!fFillErrorBand)
01964       return;
01965 
01966    // function fitting
01967    if (fFitFunctionIndexX < 0)
01968       return;
01969 
01970    // loop over all possible x values ...
01971    if (fErrorBandContinuous) {
01972       double x = 0;
01973       for (int ix = 0; ix < fErrorBandNbinsX; ix++) {
01974          // calculate x
01975          x = fErrorBandXY->GetXaxis()->GetBinCenter(ix + 1);
01976 
01977          // calculate y
01978          std::vector <double> xvec;
01979          xvec.push_back(x);
01980 
01981          // loop over all chains
01982          for (int ichain = 0; ichain < MCMCGetNChains(); ++ichain) {
01983             // calculate y
01984             double y = FitFunction(xvec, MCMCGetx(ichain));
01985 
01986             // fill histogram
01987             fErrorBandXY->Fill(x, y);
01988          }
01989 
01990          xvec.clear();
01991       }
01992    }
01993    // ... or evaluate at the data point x-values
01994    else {
01995       int ndatapoints = int(fErrorBandX.size());
01996       double x = 0;
01997 
01998       for (int ix = 0; ix < ndatapoints; ++ix) {
01999          // calculate x
02000          x = fErrorBandX.at(ix);
02001 
02002          // calculate y
02003          std::vector <double> xvec;
02004          xvec.push_back(x);
02005 
02006          // loop over all chains
02007          for (int ichain = 0; ichain < MCMCGetNChains(); ++ichain) {
02008             // calculate y
02009             double y = FitFunction(xvec, MCMCGetx(ichain));
02010 
02011             // fill histogram
02012             fErrorBandXY->Fill(x, y);
02013          }
02014 
02015          xvec.clear();
02016       }
02017    }
02018 }
02019 
02020 // *********************************************
02021 void BCIntegrate::SAInitialize()
02022 {
02023    fSAx.clear();
02024    fSAx.assign(fNvar, 0.0);
02025 }
02026 
02027 // *********************************************
02028 

Generated on Mon Aug 30 2010 22:14:54 for Bayesian Analysis Toolkit by  doxygen 1.7.1