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

BCIntegrate.h

Go to the documentation of this file.
00001 #ifndef __BCINTEGRATE__H
00002 #define __BCINTEGRATE__H
00003 
00023 // ---------------------------------------------------------
00024 
00025 #include <vector>
00026 
00027 #include <math.h>
00028 
00029 #include "BCEngineMCMC.h"
00030 #include "BCParameter.h"
00031 #include "BCDataPoint.h"
00032 
00033 // ROOT classes
00034 class TH1D;
00035 class TH2D;
00036 class TRandom3;
00037 class TMinuit;
00038 class TTree;
00039 
00040 
00041 #define BC_DEBUG 0
00042 
00043 // ---------------------------------------------------------
00044 
00045 class BCIntegrate : public BCEngineMCMC
00046 {
00047 
00048  public:
00049 
00055       enum BCIntegrationMethod { kIntMonteCarlo, kIntImportance, kIntMetropolis, kIntCuba, NIntMethods };
00056 
00059       enum BCMarginalizationMethod { kMargMonteCarlo, kMargMetropolis, NMargMethods };
00060 
00063       enum BCOptimizationMethod { kOptSA, kOptMetropolis, kOptMinuit, NOptMethods };
00064 
00067       enum BCSASchedule { kSACauchy, kSABoltzmann, kSACustom, NSAMethods };
00068 
00071       enum BCCubaMethod { kCubaVegas, kCubaSuave, kCubaDivonne, kCubaCuhre };
00072 
00080       BCIntegrate();
00081 
00084       BCIntegrate(BCParameterSet * par);
00085 
00088       BCIntegrate(const BCIntegrate & bcintegrate);
00089 
00092       virtual ~BCIntegrate();
00093 
00100       BCIntegrate & operator = (const BCIntegrate & bcintegrate);
00101 
00108       BCIntegrate::BCIntegrationMethod GetIntegrationMethod()
00109          { return fIntegrationMethod; }
00110 
00113       BCIntegrate::BCMarginalizationMethod GetMarginalizationMethod()
00114          { return fMarginalizationMethod; }
00115 
00118       BCIntegrate::BCOptimizationMethod GetOptimizationMethod()
00119          { return fOptimizationMethod; }
00120 
00123       BCIntegrate::BCOptimizationMethod GetOptimizationMethodMode()
00124          { return fOptimizationMethodMode; }
00125 
00128       BCIntegrate::BCSASchedule GetSASchedule()
00129          { return fSASchedule; }
00130 
00134       void GetRandomVector(std::vector<double> &x);
00135 
00136       virtual void GetRandomVectorMetro(std::vector<double> &x);
00137 
00143       double GetRandomPoint(std::vector<double> &x);
00144 
00150       double GetRandomPointImportance(std::vector<double> &x);
00151 
00156       void GetRandomPointMetro(std::vector<double> &x);
00157 
00162       void GetRandomPointSamplingMetro(std::vector<double> &x);
00163 
00166       int GetNiterationsPerDimension()
00167          { return fNiterPerDimension; }
00168 
00171       int GetNSamplesPer2DBin()
00172          { return fNSamplesPer2DBin; }
00173 
00176       int GetNvar()
00177          { return fNvar; }
00178 
00181       int GetNIterationsMax()
00182          { return fNIterationsMax; }
00183 
00186       int GetNIterations()
00187          { return fNIterations; }
00188 
00191       double GetRelativePrecision()
00192          { return fRelativePrecision; }
00193 
00196       double GetAbsolutePrecision()
00197          { return fAbsolutePrecision; }
00198 
00201       BCCubaMethod GetCubaIntegrationMethod()
00202          { return fCubaIntegrationMethod; }
00203 
00206       int GetCubaMinEval()
00207          { return fCubaMinEval; }
00208 
00211       int GetCubaMaxEval()
00212          { return fCubaMaxEval; }
00213 
00216       int GetCubaVerbositylevel()
00217          { return fCubaVerbosity; }
00218 
00221       int GetCubaVegasNStart()
00222          { return fCubaVegasNStart; }
00223 
00226       int GetCubaVegasNIncrease()
00227          { return fCubaVegasNIncrease; }
00228 
00231       int GetCubaSuaveNNew()
00232          { return fCubaSuaveNNew; }
00233 
00236       double GetCubaSuaveFlatness()
00237          { return fCubaSuaveFlatness; }
00238 
00241       double GetError()
00242          { return fError; }
00243 
00246       int GetNbins()
00247          { return fNbins; }
00248 
00251       TMinuit * GetMinuit();
00252 
00255       int GetMinuitErrorFlag()
00256          { return fMinuitErrorFlag; }
00257 
00260       TTree * GetMarkovChainTree()
00261          { return fMarkovChainTree; }
00262 
00265       std::vector<double> * GetMarkovChainPoint()
00266          { return &fXmetro1; }
00267 
00270       int * GetMCMCIteration()
00271          { return &fMCMCIteration; }
00272 
00275       double * GetMarkovChainValue()
00276          { return &fMarkovChainValue; }
00277 
00280       double GetSAT0()
00281          { return fSAT0; }
00282 
00285       double GetSATmin()
00286          { return fSATmin; }
00287 
00295       void SetMinuitArlist(double * arglist)
00296          { fMinuitArglist[0] = arglist[0];
00297            fMinuitArglist[1] = arglist[1]; }
00298 
00301       void SetFlagIgnorePrevOptimization(bool flag)
00302          { fFlagIgnorePrevOptimization = flag; }
00303 
00307       void SetParameters(BCParameterSet * par);
00308 
00311       void SetVarList(int * varlist);
00312 
00315       void SetVar(int index)
00316          {fVarlist[index]=1;}
00317 
00320       void SetIntegrationMethod(BCIntegrate::BCIntegrationMethod method);
00321 
00324       void SetMarginalizationMethod(BCIntegrate::BCMarginalizationMethod method)
00325          { fMarginalizationMethod = method; }
00326 
00329       void SetOptimizationMethod(BCIntegrate::BCOptimizationMethod method)
00330          { fOptimizationMethod = method; }
00331 
00335       void SetOptimizationMethodMode(BCIntegrate::BCOptimizationMethod method)
00336          { fOptimizationMethodMode = method; }
00337 
00338 
00341       void SetSASchedule(BCIntegrate::BCSASchedule schedule)
00342          { fSASchedule = schedule; }
00343 
00346       void SetNiterationsPerDimension(int niterations)
00347          { fNiterPerDimension = niterations; }
00348 
00352       void SetNSamplesPer2DBin(int n)
00353          { fNSamplesPer2DBin = n; }
00354 
00357       void SetNIterationsMax(int niterations)
00358          { fNIterationsMax = niterations; }
00359 
00363       void SetRelativePrecision(double relprecision)
00364          { fRelativePrecision = relprecision; }
00365 
00368       void SetAbsolutePrecision(double absprecision)
00369          { fAbsolutePrecision = absprecision; }
00370 
00373       void SetCubaIntegrationMethod(BCCubaMethod type);
00374 
00377       void SetCubaMinEval(int n)
00378          { fCubaMinEval = n; }
00379 
00382       void SetCubaMaxEval(int n)
00383          { fCubaMaxEval = n; }
00384 
00387       void SetCubaVerbosityLevel(int n)
00388          { fCubaVerbosity = n; }
00389 
00392       void SetCubaVegasNStart(int n)
00393          { fCubaVegasNStart = n; }
00394 
00397       void SetCubaVegasNIncrease(int n)
00398          { fCubaVegasNIncrease = n; }
00399 
00402       void SetCubaSuaveNNew(int n)
00403          { fCubaSuaveNNew = n; }
00404 
00407       void SetCubaSuaveFlatness(double p)
00408          { fCubaSuaveFlatness = p; }
00409 
00414       void SetNbins(int nbins, int index = -1);
00415 
00419       void SetFillErrorBand(bool flag = true)
00420          { fFillErrorBand=flag; }
00421 
00425       void UnsetFillErrorBand()
00426          { fFillErrorBand=false; }
00427 
00431       void SetFitFunctionIndexX(int index)
00432          { fFitFunctionIndexX = index; }
00433 
00437       void SetFitFunctionIndexY(int index)
00438          { fFitFunctionIndexY = index; }
00439 
00444       void SetFitFunctionIndices(int indexx, int indexy)
00445          { SetFitFunctionIndexX(indexx);
00446            SetFitFunctionIndexY(indexy); }
00447 
00451       void SetDataPointLowerBoundaries(BCDataPoint * datasetlowerboundaries)
00452          { fDataPointLowerBoundaries = datasetlowerboundaries; }
00453 
00457       void SetDataPointUpperBoundaries(BCDataPoint * datasetupperboundaries)
00458          { fDataPointUpperBoundaries = datasetupperboundaries; }
00459 
00463       void SetDataPointLowerBoundary(int index, double lowerboundary)
00464          { fDataPointLowerBoundaries -> SetValue(index, lowerboundary); }
00465 
00469       void SetDataPointUpperBoundary(int index, double upperboundary)
00470          { fDataPointUpperBoundaries -> SetValue(index, upperboundary); }
00471 
00474       void WriteMarkovChain(bool flag)
00475          { fFlagWriteMarkovChain = flag;
00476            fMCMCFlagWriteChainToFile = flag;
00477            fMCMCFlagWritePreRunToFile = flag; }
00478 
00481       void SetMarkovChainTree(TTree * tree)
00482          { fMarkovChainTree = tree; }
00483 
00486       void SetMarkovChainInitialPosition(std::vector<double> position);
00487 
00490       void SetMarkovChainStepSize(double stepsize)
00491          { fMarkovChainStepSize = stepsize; }
00492 
00495       void SetMarkovChainNIterations(int niterations)
00496          { fMarkovChainNIterations = niterations;
00497            fMarkovChainAutoN = false; }
00498 
00501       void SetMarkovChainAutoN(bool flag)
00502          { fMarkovChainAutoN = flag; }
00503 
00506       void SetMode(std::vector<double> mode);
00507 
00510       void SetErrorBandHisto(TH2D * h)
00511          { fErrorBandXY = h; }
00512 
00515       void SetSAT0(double T0)
00516          { fSAT0 = T0; }
00517 
00520       void SetSATmin(double Tmin)
00521          { fSATmin = Tmin; }
00522 
00523       void SetFlagWriteSAToFile(bool flag)
00524          { fFlagWriteSAToFile = flag; }
00525 
00528       void SetSATree(TTree * tree)
00529          { fTreeSA = tree; }
00530 
00533       TTree * GetSATree()
00534          { return fTreeSA; }
00535 
00538       void InitializeSATree();
00539 
00547       void DeleteVarList();
00548 
00552       void ResetVarlist(int v);
00553 
00557       void UnsetVar(int index)
00558          { fVarlist[index] = 0; }
00559 
00565       virtual double Eval(const std::vector<double> &x);
00566 
00572       virtual double LogEval(const std::vector<double> &x);
00573 
00579       virtual double EvalSampling(const std::vector<double> &x);
00580 
00586       double LogEvalSampling(const std::vector<double> &x);
00587 
00593       virtual double FitFunction(const std::vector<double> &/*x*/, const std::vector<double> &/*parameters*/)
00594          { return 0.; }
00595 
00599       double Integrate();
00600 
00606       double IntegralMC(const std::vector<double> &x, int * varlist);
00607 
00612       double IntegralMC(const std::vector<double> &x);
00613 
00618       double IntegralMetro(const std::vector<double> &x);
00619 
00624       double IntegralImportance(const std::vector<double> &x);
00625 
00632       double CubaIntegrate(BCIntegrate::BCCubaMethod method, std::vector<double> parameters_double, std::vector<double> parameters_int);
00633 
00637       double CubaIntegrate();
00638 
00646       static int CubaIntegrand(const int * ndim, const double xx[], const int * ncomp, double ff[], void *userdata);
00647 
00652       TH1D * Marginalize(BCParameter * parameter);
00653 
00659       TH2D * Marginalize(BCParameter * parameter1, BCParameter * parameter2);
00660 
00666       TH1D * MarginalizeByIntegrate(BCParameter * parameter);
00667 
00674       TH2D * MarginalizeByIntegrate(BCParameter * parameter1, BCParameter * parameter2);
00675 
00681       TH1D * MarginalizeByMetro(BCParameter * parameter);
00682 
00688       TH2D * MarginalizeByMetro(BCParameter * parameter1, BCParameter * parameter2);
00689 
00695       int MarginalizeAllByMetro(const char * name);
00696 
00700       TH1D * GetH1D(int parIndex);
00701 
00707       int GetH2DIndex(int parIndex1, int parIndex2);
00708 
00714       TH2D * GetH2D(int parIndex1, int parIndex2);
00715 
00718       void InitMetro();
00719 
00722       void SAInitialize();
00723 
00726 //      void FindMode();
00727 
00733       virtual void FindModeMinuit(std::vector<double> start = std::vector<double>(0), int printlevel = 1);
00734 
00737       void FindModeMCMC();
00738 
00744       void FindModeSA(std::vector<double> start = std::vector<double>(0));
00745 
00751       double SATemperature(double t);
00752 
00757       double SATemperatureBoltzmann(double t);
00758 
00763       double SATemperatureCauchy(double t);
00764 
00770       virtual double SATemperatureCustom(double t);
00771 
00779       std::vector<double> GetProposalPointSA(const std::vector<double> &x, int t);
00780 
00787       std::vector<double> GetProposalPointSABoltzmann(const std::vector<double> &x, int t);
00788 
00795       std::vector<double> GetProposalPointSACauchy(const std::vector<double> &x, int t);
00796 
00804       virtual std::vector<double> GetProposalPointSACustom(const std::vector<double> &x, int t);
00805 
00810       std::vector<double> SAHelperGetRandomPointOnHypersphere();
00811 
00815       double SAHelperGetRadialCauchy();
00816 
00820       double SAHelperSinusToNIntegral(int dim, double theta);
00821 
00822 
00823       static void FCNLikelihood(int &npar, double * grad, double &fval, double * par, int flag);
00824 
00828       virtual void MCMCUserIterationInterface()
00829          {}
00830 
00834       int IntegrateResetResults();
00835 
00840       std::string DumpIntegrationMethod(BCIntegrationMethod type);
00841 
00845       std::string DumpIntegrationMethod()
00846          { return DumpIntegrationMethod(fIntegrationMethod); }
00847 
00852       std::string DumpMarginalizationMethod(BCMarginalizationMethod type);
00853 
00857       std::string DumpMarginalizationMethod()
00858          { return DumpMarginalizationMethod(fMarginalizationMethod); }
00859 
00864       std::string DumpOptimizationMethod(BCOptimizationMethod type);
00865 
00869       std::string DumpOptimizationMethod()
00870          { return DumpOptimizationMethod(fOptimizationMethod); }
00871 
00875       std::string DumpUsedOptimizationMethod()
00876          { return DumpOptimizationMethod(fOptimizationMethodMode); }
00877 
00882       std::string DumpCubaIntegrationMethod(BCCubaMethod type);
00883 
00887       std::string DumpCubaIntegrationMethod()
00888          { return DumpCubaIntegrationMethod(fCubaIntegrationMethod); }
00889 
00890 
00893    protected:
00894 
00897       int fNvar;
00898 
00901       int fNbins;
00902 
00906       int fNSamplesPer2DBin;
00907 
00910       double fMarkovChainStepSize;
00911 
00912       int fMarkovChainNIterations;
00913 
00914       bool fMarkovChainAutoN;
00915 
00918       BCDataPoint * fDataPointLowerBoundaries;
00919 
00922       BCDataPoint * fDataPointUpperBoundaries;
00923 
00924       std::vector<bool> fDataFixedValues;
00925 
00929       std::vector<double> fBestFitParameters;
00930       std::vector<double> fBestFitParameterErrors;
00931 
00934       std::vector<double> fBestFitParametersMarginalized;
00935 
00938       std::vector<TH1D *> fHProb1D;
00939 
00942       std::vector<TH2D *> fHProb2D;
00943 
00946       bool fFillErrorBand;
00947 
00950       int fFitFunctionIndexX;
00951       int fFitFunctionIndexY;
00952 
00955       bool fErrorBandContinuous;
00956       std::vector<double> fErrorBandX;
00957 
00960       TH2D * fErrorBandXY;
00961 
00964       int fErrorBandNbinsX;
00965 
00968       int fErrorBandNbinsY;
00969 
00972       TMinuit * fMinuit;
00973 
00974       double fMinuitArglist[2];
00975       int fMinuitErrorFlag;
00976 
00979       double fFlagIgnorePrevOptimization;
00980 
00983       bool fFlagWriteMarkovChain;
00984 
00987       TTree * fMarkovChainTree;
00988 
00991       int fMCMCIteration;
00992 
00995       double fSAT0;
00996 
00999       double fSATmin;
01000 
01003       TTree * fTreeSA;
01004 
01007       bool fFlagWriteSAToFile;
01008 
01009       int fSANIterations;
01010       double fSATemperature;
01011       double fSALogProb;
01012       std::vector<double> fSAx;
01013 
01014    private:
01015 
01018       BCParameterSet * fx;
01019 
01022       double * fMin;
01023 
01026       double * fMax;
01027 
01030       int * fVarlist;
01031 
01034       int fNiterPerDimension;
01035 
01038       BCIntegrate::BCIntegrationMethod fIntegrationMethod;
01039 
01042       BCIntegrate::BCMarginalizationMethod fMarginalizationMethod;
01043 
01046       BCIntegrate::BCOptimizationMethod fOptimizationMethod;
01047 
01051       BCIntegrate::BCOptimizationMethod fOptimizationMethodMode;
01052 
01055       BCIntegrate::BCSASchedule fSASchedule;
01056 
01059       int fNIterationsMax;
01060 
01063       int fNIterations;
01064 
01066       double fRelativePrecision;
01067 
01069       double fAbsolutePrecision;
01070 
01072       BCCubaMethod fCubaIntegrationMethod;
01073 
01075       int fCubaMinEval;
01076 
01078       int fCubaMaxEval;
01079 
01081       int fCubaVerbosity;
01082 
01084       int fCubaVegasNStart;
01085 
01087       int fCubaVegasNIncrease;
01088 
01090       int fCubaSuaveNNew;
01091 
01093       double fCubaSuaveFlatness;
01094 
01097       double fError;
01098 
01101       int fNmetro;
01102       int fNacceptedMCMC;
01103 
01106       std::vector<double> fXmetro0;
01107 
01110       std::vector<double> fXmetro1;
01111 
01114       double fMarkovChainValue;
01115 
01118       void MCMCIterationInterface();
01119 
01122       void MCMCFillErrorBand();
01123 };
01124 
01125 // ---------------------------------------------------------
01126 
01127 #endif

Generated by  doxygen 1.7.1