00001 #ifndef __BCINTEGRATE__H
00002 #define __BCINTEGRATE__H
00003
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025 #include <vector>
00026
00027 #include <math.h>
00028
00029 #include "BAT/BCEngineMCMC.h"
00030 #include "BAT/BCParameter.h"
00031 #include "BAT/BCDataPoint.h"
00032
00033
00034 class TH1D;
00035 class TH2D;
00036 class TRandom3;
00037 class TMinuit;
00038 class TTree;
00039
00040
00041 #define DEBUG 0
00042
00043
00044
00045 class BCIntegrate : public BCEngineMCMC
00046 {
00047
00048 public:
00049
00051
00052
00055 enum BCIntegrationMethod { kIntMonteCarlo, kIntImportance, kIntMetropolis, kIntCuba };
00056
00059 enum BCMarginalizationMethod { kMargMonteCarlo, kMargMetropolis };
00060
00063 enum BCOptimizationMethod { kOptSimulatedAnnealing, kOptMetropolis, kOptMinuit };
00064
00065
00066
00068
00069
00072 BCIntegrate();
00073
00076 BCIntegrate(BCParameterSet * par);
00077
00080 virtual ~BCIntegrate();
00081
00082
00083
00085
00086
00089 BCIntegrate::BCIntegrationMethod GetIntegrationMethod()
00090 { return fIntegrationMethod; };
00091
00094 BCIntegrate::BCMarginalizationMethod GetMarginalizationMethod()
00095 { return fMarginalizationMethod; };
00096
00099 BCIntegrate::BCOptimizationMethod GetOptimizationMethod()
00100 { return fOptimizationMethod; };
00101
00105 void GetRandomVector(std::vector <double> &x);
00106
00107 virtual void GetRandomVectorMetro(std::vector <double> &x);
00108
00114 double GetRandomPoint(std::vector <double> &x);
00115
00121 double GetRandomPointImportance(std::vector <double> &x);
00122
00127 void GetRandomPointMetro(std::vector <double> &x);
00128
00133 void GetRandomPointSamplingMetro(std::vector <double> &x);
00134
00137 int GetNiterationsPerDimension()
00138 { return fNiterPerDimension; };
00139
00142 int GetNSamplesPer2DBin()
00143 { return fNSamplesPer2DBin; };
00144
00147 int GetNvar()
00148 { return fNvar; };
00149
00152 int GetNIterationsMax()
00153 { return fNIterationsMax; };
00154
00157 int GetNIterations()
00158 { return fNIterations; };
00159
00162 double GetRelativePrecision()
00163 { return fRelativePrecision; };
00164
00165
00166
00167 double GetError()
00168 { return fError; };
00169
00170
00171
00172 int GetNbins()
00173 { return fNbins; };
00174
00175
00176
00177 TMinuit * GetMinuit();
00178
00179 int GetMinuitErrorFlag()
00180 { return fMinuitErrorFlag; };
00181
00182
00183
00184 TTree * GetMarkovChainTree()
00185 { return fMarkovChainTree; };
00186
00187
00188
00189 std::vector<double> * GetMarkovChainPoint()
00190 { return &fXmetro1; };
00191
00194 int * GetMCMCIteration()
00195 { return &fMCMCIteration; };
00196
00199 double * GetMarkovChainValue()
00200 { return &fMarkovChainValue; };
00201
00202
00203
00205
00206
00207 void SetMinuitArlist(double * arglist)
00208 { fMinuitArglist[0] = arglist[0];
00209 fMinuitArglist[1] = arglist[1]; };
00210
00214 void SetParameters(BCParameterSet * par);
00215
00218 void SetVarList(int * varlist);
00219
00222 void SetVar(int index){fVarlist[index]=1;};
00223
00226 void SetIntegrationMethod(BCIntegrate::BCIntegrationMethod method);
00227
00230 void SetMarginalizationMethod(BCIntegrate::BCMarginalizationMethod method)
00231 { fMarginalizationMethod = method; };
00232
00235 void SetOptimizationMethod(BCIntegrate::BCOptimizationMethod method)
00236 { fOptimizationMethod = method; };
00237
00240 void SetNiterationsPerDimension(int niterations)
00241 { fNiterPerDimension = niterations; };
00242
00246 void SetNSamplesPer2DBin(int n)
00247 { fNSamplesPer2DBin = n; };
00248
00251 void SetNIterationsMax(int niterations)
00252 { fNIterationsMax = niterations; };
00253
00257 void SetRelativePrecision(double relprecision)
00258 { fRelativePrecision = relprecision; };
00259
00269 void SetNbins(int nbins, int index = -1);
00270
00271
00272
00273
00274
00275
00276
00277
00278 void SetFillErrorBand(bool flag = true)
00279 { fFillErrorBand=flag; };
00280
00281
00282
00283
00284 void UnsetFillErrorBand()
00285 { fFillErrorBand=false; };
00286
00290 void SetFitFunctionIndexX(int index)
00291 { fFitFunctionIndexX = index; };
00292
00296 void SetFitFunctionIndexY(int index)
00297 { fFitFunctionIndexY = index; };
00298
00299 void SetFitFunctionIndices(int indexx, int indexy)
00300 { this -> SetFitFunctionIndexX(indexx);
00301 this -> SetFitFunctionIndexY(indexy); };
00302
00306 void SetDataPointLowerBoundaries(BCDataPoint * datasetlowerboundaries)
00307 { fDataPointLowerBoundaries = datasetlowerboundaries; };
00308
00312 void SetDataPointUpperBoundaries(BCDataPoint * datasetupperboundaries)
00313 { fDataPointUpperBoundaries = datasetupperboundaries; };
00314
00318 void SetDataPointLowerBoundary(int index, double lowerboundary)
00319 { fDataPointLowerBoundaries -> SetValue(index, lowerboundary); };
00320
00324 void SetDataPointUpperBoundary(int index, double upperboundary)
00325 { fDataPointUpperBoundaries -> SetValue(index, upperboundary); };
00326
00330 void WriteMarkovChain(bool flag)
00331 { fFlagWriteMarkovChain = flag; fMCMCFlagWriteChainToFile = flag;};
00332
00335 void SetMarkovChainTree(TTree * tree)
00336 { fMarkovChainTree = tree; };
00337
00340 void SetMarkovChainInitialPosition(std::vector<double> position);
00341
00344 void SetMarkovChainStepSize(double stepsize)
00345 { fMarkovChainStepSize = stepsize; };
00346
00349 void SetMarkovChainNIterations(int niterations)
00350 { fMarkovChainNIterations = niterations;
00351 fMarkovChainAutoN = false; }
00352
00355 void SetMarkovChainAutoN(bool flag)
00356 { fMarkovChainAutoN = flag; };
00357
00360 void SetMode(std::vector <double> mode);
00361
00364 void SetErrorBandHisto(TH2D * h)
00365 { fErrorBandXY = h; };
00366
00367
00368
00370
00371
00374 void DeleteVarList();
00375
00379 void ResetVarlist(int v);
00380
00384 void UnsetVar(int index)
00385 { fVarlist[index] = 0; };
00386
00392 virtual double Eval(std::vector <double> x);
00393
00399 virtual double LogEval(std::vector <double> x);
00400
00406 virtual double EvalSampling(std::vector <double> x);
00407
00413 double LogEvalSampling(std::vector <double> x);
00414
00421 double EvalPrint(std::vector <double> x);
00422
00428 virtual double FitFunction(std::vector <double> x, std::vector <double> parameters)
00429 { return 0.; };
00430
00434 double Integrate();
00435
00441 double IntegralMC(std::vector <double> x, int * varlist);
00442
00443 double IntegralMC(std::vector <double> x);
00444
00449 double IntegralMetro(std::vector <double> x);
00450
00455 double IntegralImportance(std::vector <double> x);
00456
00463 double CubaIntegrate(int method, std::vector<double> parameters_double, std::vector<double> parameters_int);
00464
00465 double CubaIntegrate();
00466
00474 static void CubaIntegrand(const int * ndim, const double xx[], const int * ncomp, double ff[]);
00475
00480 TH1D * Marginalize(BCParameter * parameter);
00481
00487 TH2D * Marginalize(BCParameter * parameter1, BCParameter * parameter2);
00488
00494 TH1D * MarginalizeByIntegrate(BCParameter * parameter);
00495
00502 TH2D * MarginalizeByIntegrate(BCParameter * parameter1, BCParameter * parameter2);
00503
00509 TH1D * MarginalizeByMetro(BCParameter * parameter);
00510
00516 TH2D * MarginalizeByMetro(BCParameter * parameter1, BCParameter * parameter2);
00517
00523 int MarginalizeAllByMetro(const char * name);
00524
00529 TH1D * GetH1D(int parIndex);
00530
00536 int GetH2DIndex(int parIndex1, int parIndex2);
00537
00543 TH2D * GetH2D(int parIndex1, int parIndex2);
00544
00547 void InitMetro();
00548
00551
00552
00558 void FindModeMinuit(std::vector<double> start = std::vector<double>(0), int printlevel = 1);
00559
00562 void FindModeMCMC();
00563
00564
00565 static void FCNLikelihood(int &npar, double * grad, double &fval, double * par, int flag);
00566
00567
00568
00569
00570 virtual void MCMCUserIterationInterface()
00571 {};
00572
00573
00574
00575 private:
00576
00579 BCParameterSet * fx;
00580
00583 double * fMin;
00584
00587 double * fMax;
00588
00591 int * fVarlist;
00592
00595 int fNiterPerDimension;
00596
00599 BCIntegrate::BCIntegrationMethod fIntegrationMethod;
00600
00603 BCIntegrate::BCMarginalizationMethod fMarginalizationMethod;
00604
00607 BCIntegrate::BCOptimizationMethod fOptimizationMethod;
00608
00611 int fNIterationsMax;
00612
00615 int fNIterations;
00616
00619 double fRelativePrecision;
00620
00623 double fError;
00624
00627 int fNmetro;
00628 int fNacceptedMCMC;
00629
00632 std::vector <double> fXmetro0;
00633
00636 std::vector <double> fXmetro1;
00637
00640 double fMarkovChainValue;
00641
00642
00643
00644 void MCMCIterationInterface();
00645
00646
00647
00648 void MCMCFillErrorBand();
00649
00650 protected:
00651
00654 int fNvar;
00655
00658 int fNbins;
00659
00663 int fNSamplesPer2DBin;
00664
00667 double fMarkovChainStepSize;
00668
00669 int fMarkovChainNIterations;
00670
00671 bool fMarkovChainAutoN;
00672
00675 BCDataPoint * fDataPointLowerBoundaries;
00676
00679 BCDataPoint * fDataPointUpperBoundaries;
00680
00681 std::vector <bool> fDataFixedValues;
00682
00685 TRandom3 * fRandom;
00686
00689 std::vector <double> fBestFitParameters;
00690
00693 std::vector <double> fBestFitParametersMarginalized;
00694
00697 std::vector <TH1D *> fHProb1D;
00698
00701 std::vector <TH2D *> fHProb2D;
00702
00705 bool fFillErrorBand;
00706
00709 int fFitFunctionIndexX;
00710 int fFitFunctionIndexY;
00711
00714 bool fErrorBandContinuous;
00715 std::vector <double> fErrorBandX;
00716
00719 TH2D * fErrorBandXY;
00720
00723 int fErrorBandNbinsX;
00724
00727 int fErrorBandNbinsY;
00728
00731 TMinuit * fMinuit;
00732
00733 double fMinuitArglist[2];
00734 int fMinuitErrorFlag;
00735
00738 bool fFlagWriteMarkovChain;
00739
00742 TTree * fMarkovChainTree;
00743
00746 int fMCMCIteration;
00747
00748 };
00749
00750
00751
00752 #endif