00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "BCEngineMCMC.h"
00011
00012 #include "BCLog.h"
00013 #include "BCMath.h"
00014 #include "BCParameter.h"
00015
00016 #include <TH1D.h>
00017 #include <TH2D.h>
00018 #include <TTree.h>
00019 #include <TRandom3.h>
00020
00021 #include <math.h>
00022
00023
00024 BCEngineMCMC::BCEngineMCMC()
00025 {
00026
00027 MCMCSetValuesDefault();
00028
00029
00030 fRandom = new TRandom3(0);
00031 }
00032
00033
00034 BCEngineMCMC::BCEngineMCMC(int n)
00035 {
00036
00037 fMCMCNChains = n;
00038
00039
00040 BCEngineMCMC();
00041 }
00042
00043
00044 void BCEngineMCMC::MCMCSetValuesDefault()
00045 {
00046 fMCMCNParameters = 0;
00047 fMCMCFlagWriteChainToFile = false;
00048 fMCMCFlagWritePreRunToFile = false;
00049 fMCMCFlagPreRun = true;
00050 fMCMCFlagRun = false;
00051 fMCMCFlagFillHistograms = true;
00052 fMCMCEfficiencyMin = 0.15;
00053 fMCMCEfficiencyMax = 0.50;
00054 fMCMCFlagInitialPosition = 1;
00055 fMCMCNLag = 1;
00056 fMCMCCurrentIteration = -1;
00057 fMCMCCurrentChain = -1;
00058
00059 MCMCSetValuesDetail();
00060 }
00061
00062
00063 void BCEngineMCMC::MCMCSetValuesQuick()
00064 {
00065 fMCMCNChains = 2;
00066 fMCMCNIterationsMax = 10000;
00067 fMCMCNIterationsRun = 10000;
00068 fMCMCNIterationsPreRunMin = 500;
00069 fMCMCFlagInitialPosition = 1;
00070 fMCMCRValueUseStrict = false;
00071 fMCMCRValueCriterion = 0.1;
00072 fMCMCRValueParametersCriterion = 0.1;
00073 fMCMCNIterationsConvergenceGlobal = -1;
00074 fMCMCFlagConvergenceGlobal = false;
00075 fMCMCRValue = 100;
00076 fMCMCNIterationsUpdate = 1000;
00077 fMCMCNIterationsUpdateMax = 10000;
00078 fMCMCFlagOrderParameters = true;
00079 fMCMCCurrentIteration = -1;
00080 fMCMCCurrentChain = -1;
00081 }
00082
00083
00084 void BCEngineMCMC::MCMCSetValuesDetail()
00085 {
00086 fMCMCNChains = 5;
00087 fMCMCNIterationsMax = 1000000;
00088 fMCMCNIterationsRun = 100000;
00089 fMCMCNIterationsPreRunMin = 500;
00090 fMCMCFlagInitialPosition = 1;
00091 fMCMCRValueUseStrict = false;
00092 fMCMCRValueCriterion = 0.1;
00093 fMCMCRValueParametersCriterion = 0.1;
00094 fMCMCNIterationsConvergenceGlobal = -1;
00095 fMCMCFlagConvergenceGlobal = false;
00096 fMCMCRValue = 100;
00097 fMCMCNIterationsUpdate = 1000;
00098 fMCMCNIterationsUpdateMax = 10000;
00099 fMCMCFlagOrderParameters = true;
00100 fMCMCCurrentIteration = -1;
00101 fMCMCCurrentChain = -1;
00102 }
00103
00104
00105 void BCEngineMCMC::MCMCSetPrecision(BCEngineMCMC::Precision precision)
00106 {
00107 switch(precision) {
00108 case BCEngineMCMC::kLow:
00109 fMCMCNChains = 1;
00110 fMCMCNLag = 1;
00111 fMCMCNIterationsMax = 10000;
00112 fMCMCNIterationsRun = 10000;
00113 fMCMCNIterationsPreRunMin = 100;
00114 fMCMCNIterationsUpdate = 1000;
00115 fMCMCNIterationsUpdateMax = 10000;
00116 fMCMCRValueCriterion = 0.1;
00117 fMCMCRValueParametersCriterion = 0.1;
00118 fMCMCRValue = 100;
00119 break;
00120 case BCEngineMCMC::kMedium:
00121 fMCMCNChains = 5;
00122 fMCMCNLag = 1;
00123 fMCMCNIterationsMax = 100000;
00124 fMCMCNIterationsRun = 100000;
00125 fMCMCNIterationsPreRunMin = 100;
00126 fMCMCNIterationsUpdate = 1000;
00127 fMCMCNIterationsUpdateMax = 10000;
00128 fMCMCRValueCriterion = 0.1;
00129 fMCMCRValueParametersCriterion = 0.1;
00130 fMCMCRValue = 100;
00131 break;
00132 case BCEngineMCMC::kHigh:
00133 fMCMCNChains = 10;
00134 fMCMCNLag = 10;
00135 fMCMCNIterationsMax = 1000000;
00136 fMCMCNIterationsRun = 1000000;
00137 fMCMCNIterationsPreRunMin = 100;
00138 fMCMCNIterationsUpdate = 1000;
00139 fMCMCNIterationsUpdateMax = 10000;
00140 fMCMCRValueCriterion = 0.1;
00141 fMCMCRValueParametersCriterion = 0.1;
00142 fMCMCRValue = 100;
00143 break;
00144 case BCEngineMCMC::kVeryHigh:
00145 fMCMCNChains = 10;
00146 fMCMCNLag = 10;
00147 fMCMCNIterationsMax = 10000000;
00148 fMCMCNIterationsRun = 10000000;
00149 fMCMCNIterationsPreRunMin = 100;
00150 fMCMCNIterationsUpdate = 1000;
00151 fMCMCNIterationsUpdateMax = 10000;
00152 fMCMCRValueCriterion = 0.1;
00153 fMCMCRValueParametersCriterion = 0.1;
00154 fMCMCRValue = 100;
00155 break;
00156 }
00157
00158
00159 MCMCInitialize();
00160 }
00161
00162
00163 BCEngineMCMC::~BCEngineMCMC()
00164 {
00165
00166 if (fRandom)
00167 delete fRandom;
00168
00169
00170 for (int i = 0; i < int(fMCMCH1Marginalized.size()); ++i)
00171 if (fMCMCH1Marginalized[i])
00172 delete fMCMCH1Marginalized[i];
00173 fMCMCH1Marginalized.clear();
00174
00175
00176 for (int i = 0; i < int(fMCMCH2Marginalized.size()); ++i)
00177 if (fMCMCH2Marginalized[i])
00178 delete fMCMCH2Marginalized[i];
00179 fMCMCH2Marginalized.clear();
00180 }
00181
00182
00183 BCEngineMCMC::BCEngineMCMC(const BCEngineMCMC & enginemcmc)
00184 {
00185 fMCMCPointerToGetProposalPoint = enginemcmc.fMCMCPointerToGetProposalPoint;
00186 fMCMCNParameters = enginemcmc.fMCMCNParameters;
00187 fMCMCBoundaryMin = enginemcmc.fMCMCBoundaryMin;
00188 fMCMCBoundaryMax = enginemcmc.fMCMCBoundaryMax;
00189 fMCMCFlagsFillHistograms = enginemcmc.fMCMCFlagsFillHistograms;
00190 fMCMCNChains = enginemcmc.fMCMCNChains;
00191 fMCMCNLag = enginemcmc.fMCMCNLag;
00192 fMCMCNIterations = enginemcmc.fMCMCNIterations;
00193 fMCMCCurrentIteration = enginemcmc.fMCMCCurrentIteration;
00194 fMCMCCurrentChain = enginemcmc.fMCMCCurrentChain;
00195 fMCMCNIterationsUpdate = enginemcmc.fMCMCNIterationsUpdate;
00196 fMCMCNIterationsUpdateMax = enginemcmc.fMCMCNIterationsUpdateMax;
00197 fMCMCNIterationsConvergenceGlobal = enginemcmc.fMCMCNIterationsConvergenceGlobal;
00198 fMCMCFlagConvergenceGlobal = enginemcmc.fMCMCFlagConvergenceGlobal;
00199 fMCMCNIterationsMax = enginemcmc.fMCMCNIterationsMax;
00200 fMCMCNIterationsRun = enginemcmc.fMCMCNIterationsRun;
00201 fMCMCNIterationsPreRunMin = enginemcmc.fMCMCNIterationsPreRunMin;
00202 fMCMCNTrialsTrue = enginemcmc.fMCMCNTrialsTrue;
00203 fMCMCNTrialsFalse = enginemcmc.fMCMCNTrialsFalse;
00204 fMCMCFlagWriteChainToFile = enginemcmc.fMCMCFlagWriteChainToFile;
00205 fMCMCFlagWritePreRunToFile = enginemcmc.fMCMCFlagWritePreRunToFile;
00206 fMCMCTrialFunctionScaleFactor = enginemcmc.fMCMCTrialFunctionScaleFactor;
00207 fMCMCTrialFunctionScaleFactorStart = enginemcmc.fMCMCTrialFunctionScaleFactorStart;
00208 fMCMCFlagPreRun = enginemcmc.fMCMCFlagPreRun;
00209 fMCMCFlagRun = enginemcmc.fMCMCFlagRun;
00210 fMCMCInitialPosition = enginemcmc.fMCMCInitialPosition;
00211 fMCMCEfficiencyMin = enginemcmc.fMCMCEfficiencyMin;
00212 fMCMCEfficiencyMax = enginemcmc.fMCMCEfficiencyMax;
00213 fMCMCFlagInitialPosition = enginemcmc.fMCMCFlagInitialPosition;
00214 fMCMCFlagOrderParameters = enginemcmc.fMCMCFlagOrderParameters;
00215 fMCMCFlagFillHistograms = enginemcmc.fMCMCFlagFillHistograms;
00216 fMCMCPhase = enginemcmc.fMCMCPhase;
00217 fMCMCCycle = enginemcmc.fMCMCCycle;
00218 fMCMCx = enginemcmc.fMCMCx;
00219 fMCMCxMax = enginemcmc.fMCMCxMax;
00220 fMCMCxMean = enginemcmc.fMCMCxMean;
00221 fMCMCxVar = enginemcmc.fMCMCxVar;
00222 fMCMCxLocal = enginemcmc.fMCMCxLocal;
00223 fMCMCprob = enginemcmc.fMCMCprob;
00224 fMCMCprobMax = enginemcmc.fMCMCprobMax;
00225 fMCMCprobMean = enginemcmc.fMCMCprobMean;
00226 fMCMCprobVar = enginemcmc.fMCMCprobVar;
00227 fMCMCRValueUseStrict = enginemcmc.fMCMCRValueUseStrict;
00228 fMCMCRValueCriterion = enginemcmc.fMCMCRValueCriterion ;
00229 fMCMCRValueParametersCriterion = enginemcmc.fMCMCRValueParametersCriterion;
00230 fMCMCRValue = enginemcmc.fMCMCRValue;
00231 fMCMCRValueParameters = enginemcmc.fMCMCRValueParameters;
00232 if (enginemcmc.fRandom)
00233 fRandom = new TRandom3(*enginemcmc.fRandom);
00234 else
00235 fRandom = 0;
00236 fMCMCH1NBins = enginemcmc.fMCMCH1NBins;
00237
00238 for (int i = 0; i < int(enginemcmc.fMCMCH1Marginalized.size()); ++i) {
00239 if (enginemcmc.fMCMCH1Marginalized.at(i))
00240 fMCMCH1Marginalized.push_back(new TH1D(*(enginemcmc.fMCMCH1Marginalized.at(i))));
00241 else
00242 fMCMCH1Marginalized.push_back(0);
00243 }
00244
00245 for (int i = 0; i < int(enginemcmc.fMCMCH2Marginalized.size()); ++i) {
00246 if (enginemcmc.fMCMCH2Marginalized.at(i))
00247 fMCMCH2Marginalized.push_back(new TH2D(*(enginemcmc.fMCMCH2Marginalized.at(i))));
00248 else
00249 fMCMCH2Marginalized.push_back(0);
00250 }
00251
00252 for (int i = 0; i < int(enginemcmc.fMCMCTrees.size()); ++i) {
00253
00254 fMCMCTrees.push_back(0);
00255 }
00256 }
00257
00258
00259 BCEngineMCMC & BCEngineMCMC::operator = (const BCEngineMCMC & enginemcmc)
00260 {
00261 fMCMCPointerToGetProposalPoint = enginemcmc.fMCMCPointerToGetProposalPoint;
00262 fMCMCNParameters = enginemcmc.fMCMCNParameters;
00263 fMCMCBoundaryMin = enginemcmc.fMCMCBoundaryMin;
00264 fMCMCBoundaryMax = enginemcmc.fMCMCBoundaryMax;
00265 fMCMCFlagsFillHistograms = enginemcmc.fMCMCFlagsFillHistograms;
00266 fMCMCNChains = enginemcmc.fMCMCNChains;
00267 fMCMCNLag = enginemcmc.fMCMCNLag;
00268 fMCMCNIterations = enginemcmc.fMCMCNIterations;
00269 fMCMCCurrentIteration = enginemcmc.fMCMCCurrentIteration;
00270 fMCMCCurrentChain = enginemcmc.fMCMCCurrentChain;
00271 fMCMCNIterationsUpdate = enginemcmc.fMCMCNIterationsUpdate;
00272 fMCMCNIterationsUpdateMax = enginemcmc.fMCMCNIterationsUpdateMax;
00273 fMCMCNIterationsConvergenceGlobal = enginemcmc.fMCMCNIterationsConvergenceGlobal;
00274 fMCMCFlagConvergenceGlobal = enginemcmc.fMCMCFlagConvergenceGlobal;
00275 fMCMCNIterationsMax = enginemcmc.fMCMCNIterationsMax;
00276 fMCMCNIterationsRun = enginemcmc.fMCMCNIterationsRun;
00277 fMCMCNIterationsPreRunMin = enginemcmc.fMCMCNIterationsPreRunMin;
00278 fMCMCNTrialsTrue = enginemcmc.fMCMCNTrialsTrue;
00279 fMCMCNTrialsFalse = enginemcmc.fMCMCNTrialsFalse;
00280 fMCMCFlagWriteChainToFile = enginemcmc.fMCMCFlagWriteChainToFile;
00281 fMCMCFlagWritePreRunToFile = enginemcmc.fMCMCFlagWritePreRunToFile;
00282 fMCMCTrialFunctionScaleFactor = enginemcmc.fMCMCTrialFunctionScaleFactor;
00283 fMCMCTrialFunctionScaleFactorStart = enginemcmc.fMCMCTrialFunctionScaleFactorStart;
00284 fMCMCFlagPreRun = enginemcmc.fMCMCFlagPreRun;
00285 fMCMCFlagRun = enginemcmc.fMCMCFlagRun;
00286 fMCMCInitialPosition = enginemcmc.fMCMCInitialPosition;
00287 fMCMCEfficiencyMin = enginemcmc.fMCMCEfficiencyMin;
00288 fMCMCEfficiencyMax = enginemcmc.fMCMCEfficiencyMax;
00289 fMCMCFlagInitialPosition = enginemcmc.fMCMCFlagInitialPosition;
00290 fMCMCFlagOrderParameters = enginemcmc.fMCMCFlagOrderParameters;
00291 fMCMCFlagFillHistograms = enginemcmc.fMCMCFlagFillHistograms;
00292 fMCMCPhase = enginemcmc.fMCMCPhase;
00293 fMCMCCycle = enginemcmc.fMCMCCycle;
00294 fMCMCx = enginemcmc.fMCMCx;
00295 fMCMCxMax = enginemcmc.fMCMCxMax;
00296 fMCMCxMean = enginemcmc.fMCMCxMean;
00297 fMCMCxVar = enginemcmc.fMCMCxVar;
00298 fMCMCxLocal = enginemcmc.fMCMCxLocal;
00299 fMCMCprob = enginemcmc.fMCMCprob;
00300 fMCMCprobMax = enginemcmc.fMCMCprobMax;
00301 fMCMCprobMean = enginemcmc.fMCMCprobMean;
00302 fMCMCprobVar = enginemcmc.fMCMCprobVar;
00303 fMCMCRValueUseStrict = enginemcmc.fMCMCRValueUseStrict;
00304 fMCMCRValueCriterion = enginemcmc.fMCMCRValueCriterion ;
00305 fMCMCRValueParametersCriterion = enginemcmc.fMCMCRValueParametersCriterion;
00306 fMCMCRValue = enginemcmc.fMCMCRValue;
00307 fMCMCRValueParameters = enginemcmc.fMCMCRValueParameters;
00308 if (enginemcmc.fRandom)
00309 fRandom = new TRandom3(*enginemcmc.fRandom);
00310 else
00311 fRandom = 0;
00312 fMCMCH1NBins = enginemcmc.fMCMCH1NBins;
00313
00314 for (int i = 0; i < int(enginemcmc.fMCMCH1Marginalized.size()); ++i) {
00315 if (enginemcmc.fMCMCH1Marginalized.at(i))
00316 fMCMCH1Marginalized.push_back(new TH1D(*(enginemcmc.fMCMCH1Marginalized.at(i))));
00317 else
00318 fMCMCH1Marginalized.push_back(0);
00319 }
00320
00321 for (int i = 0; i < int(enginemcmc.fMCMCH2Marginalized.size()); ++i) {
00322 if (enginemcmc.fMCMCH2Marginalized.at(i))
00323 fMCMCH2Marginalized.push_back(new TH2D(*(enginemcmc.fMCMCH2Marginalized.at(i))));
00324 else
00325 fMCMCH2Marginalized.push_back(0);
00326 }
00327
00328 for (int i = 0; i < int(enginemcmc.fMCMCTrees.size()); ++i) {
00329
00330 fMCMCTrees.push_back(0);
00331 }
00332
00333
00334 return *this;
00335 }
00336
00337
00338
00339 TH1D * BCEngineMCMC::MCMCGetH1Marginalized(int index)
00340 {
00341
00342 if (index < 0 || index >= int(fMCMCH1Marginalized.size()))
00343 {
00344 BCLog::OutWarning("BCEngineMCMC::MCMCGetH1Marginalized. Index out of range.");
00345 return 0;
00346 }
00347
00348 return fMCMCH1Marginalized[index];
00349 }
00350
00351
00352 TH2D * BCEngineMCMC::MCMCGetH2Marginalized(int index1, int index2)
00353 {
00354 int counter = 0;
00355 int index = 0;
00356
00357
00358 for(int i = 0; i < fMCMCNParameters; i++)
00359 for (int j = 0; j < i; ++j)
00360 {
00361 if(j == index1 && i == index2)
00362 index = counter;
00363 counter++;
00364 }
00365
00366
00367 if (index < 0 || index >= int(fMCMCH2Marginalized.size()))
00368 {
00369 BCLog::OutWarning("BCEngineMCMC::MCMCGetH2Marginalized. Index out of range.");
00370 return 0;
00371 }
00372
00373 return fMCMCH2Marginalized[index];
00374 }
00375
00376
00377 std::vector<double> BCEngineMCMC::MCMCGetMaximumPoint(int i)
00378 {
00379
00380 std::vector<double> x;
00381
00382
00383 if (i < 0 || i >= fMCMCNChains)
00384 return x;
00385
00386
00387 for (int j = 0; j < fMCMCNParameters; ++j)
00388 x.push_back(fMCMCxMax.at(i * fMCMCNParameters + j));
00389
00390 return x;
00391 }
00392
00393
00394 void BCEngineMCMC::MCMCSetNChains(int n)
00395 {
00396 fMCMCNChains = n;
00397
00398
00399 MCMCInitialize();
00400 }
00401
00402
00403 void BCEngineMCMC::MCMCSetInitialPositions(const std::vector<double> & x0s)
00404 {
00405
00406 fMCMCInitialPosition.clear();
00407
00408
00409 int n = int(x0s.size());
00410
00411 for (int i = 0; i < n; ++i)
00412 fMCMCInitialPosition.push_back(x0s.at(i));
00413
00414
00415 MCMCSetFlagInitialPosition(2);
00416 }
00417
00418
00419 void BCEngineMCMC::MCMCSetInitialPositions(std::vector< std::vector<double> > x0s)
00420 {
00421
00422 std::vector<double> y0s;
00423
00424
00425 for (int i = 0; i < int(x0s.size()); ++i)
00426 for (int j = 0; j < int((x0s.at(i)).size()); ++j)
00427 y0s.push_back((x0s.at(i)).at(j));
00428
00429 MCMCSetInitialPositions(y0s);
00430 }
00431
00432
00433 void BCEngineMCMC::MCMCSetFlagFillHistograms(bool flag)
00434 {
00435 fMCMCFlagFillHistograms = flag;
00436
00437 for (int i = 0; i < fMCMCNParameters; ++i)
00438 fMCMCFlagsFillHistograms[i] = flag;
00439 }
00440
00441
00442 void BCEngineMCMC::MCMCSetFlagFillHistograms(int index, bool flag)
00443 {
00444
00445 if (index < 0 || index > fMCMCNParameters)
00446 {
00447 BCLog::OutWarning("BCEngineMCMC :MCMCSetFlagFillHistograms. Index out of range.");
00448 return;
00449 }
00450
00451
00452 fMCMCFlagsFillHistograms[index] = flag;
00453 }
00454
00455
00456 void BCEngineMCMC::MCMCSetMarkovChainTrees(std::vector<TTree *> trees)
00457 {
00458
00459 fMCMCTrees.clear();
00460
00461
00462 for (int i = 0; i < int(trees.size()); ++i)
00463 fMCMCTrees.push_back(trees[i]);
00464 }
00465
00466
00467 void BCEngineMCMC::MCMCInitializeMarkovChainTrees()
00468 {
00469
00470 fMCMCTrees.clear();
00471
00472
00473 for (int i = 0; i < fMCMCNChains; ++i) {
00474 fMCMCTrees.push_back(new TTree(TString::Format("MarkovChainTree_%i", i), "MarkovChainTree"));
00475 fMCMCTrees[i]->Branch("Iteration", &fMCMCNIterations[i], "iteration/I");
00476 fMCMCTrees[i]->Branch("NParameters", &fMCMCNParameters, "parameters/I");
00477 fMCMCTrees[i]->Branch("LogProbability", &fMCMCprob[i], "log(probability)/D");
00478 fMCMCTrees[i]->Branch("Phase", &fMCMCPhase, "phase/I");
00479 fMCMCTrees[i]->Branch("Cycle", &fMCMCCycle, "cycle/I");
00480
00481 for (int j = 0; j < fMCMCNParameters; ++j)
00482 fMCMCTrees[i]->Branch(TString::Format("Parameter%i", j),
00483 &fMCMCx[i * fMCMCNParameters + j],
00484 TString::Format("parameter %i/D", j));
00485 }
00486 }
00487
00488
00489 void BCEngineMCMC::MCMCTrialFunction(int ichain, std::vector<double> &x)
00490 {
00491
00492 for (int i = 0; i < fMCMCNParameters; ++i)
00493 x[i] = MCMCTrialFunctionSingle(ichain, i);
00494 }
00495
00496
00497 double BCEngineMCMC::MCMCTrialFunctionSingle(int ichain, int iparameter)
00498 {
00499
00500
00501
00502
00503
00504
00505 return fRandom->BreitWigner(0.0, fMCMCTrialFunctionScaleFactor[ichain * fMCMCNParameters + iparameter]);
00506 }
00507
00508
00509 std::vector<double> BCEngineMCMC::MCMCGetTrialFunctionScaleFactor(int ichain)
00510 {
00511
00512 std::vector<double> x;
00513
00514
00515 if (ichain < 0 || ichain >= fMCMCNChains)
00516 return x;
00517
00518
00519 for (int j = 0; j < fMCMCNParameters; ++j)
00520 x.push_back(fMCMCTrialFunctionScaleFactor.at(ichain * fMCMCNParameters + j));
00521
00522 return x;
00523 }
00524
00525
00526 double BCEngineMCMC::MCMCGetTrialFunctionScaleFactor(int ichain, int ipar)
00527 {
00528
00529 if (ichain < 0 || ichain >= fMCMCNChains)
00530 return 0;
00531
00532
00533 if (ipar < 0 || ipar >= fMCMCNParameters)
00534 return 0;
00535
00536
00537 return fMCMCTrialFunctionScaleFactor.at(ichain * fMCMCNChains + ipar);
00538 }
00539
00540
00541 std::vector<double> BCEngineMCMC::MCMCGetx(int ichain)
00542 {
00543
00544 std::vector<double> x;
00545
00546
00547 if (ichain < 0 || ichain >= fMCMCNChains)
00548 return x;
00549
00550
00551 for (int j = 0; j < fMCMCNParameters; ++j)
00552 x.push_back(fMCMCx.at(ichain * fMCMCNParameters + j));
00553
00554 return x;
00555 }
00556
00557
00558 double BCEngineMCMC::MCMCGetx(int ichain, int ipar)
00559 {
00560
00561 if (ichain < 0 || ichain >= fMCMCNChains)
00562 return 0;
00563
00564
00565 if (ipar < 0 || ipar >= fMCMCNParameters)
00566 return 0;
00567
00568
00569 return fMCMCx.at(ichain * fMCMCNParameters + ipar);
00570 }
00571
00572
00573 double BCEngineMCMC::MCMCGetLogProbx(int ichain)
00574 {
00575
00576 if (ichain < 0 || ichain >= fMCMCNChains)
00577 return -1;
00578
00579
00580 return fMCMCprob.at(ichain);
00581 }
00582
00583
00584 bool BCEngineMCMC::MCMCGetProposalPointMetropolis(int chain, std::vector<double> &x)
00585 {
00586
00587 MCMCTrialFunction(chain, x);
00588
00589
00590 for (int i = 0; i < fMCMCNParameters; ++i)
00591 x[i] = fMCMCx[chain * fMCMCNParameters + i] + x[i] * (fMCMCBoundaryMax.at(i) - fMCMCBoundaryMin.at(i));
00592
00593
00594 for (int i = 0; i < fMCMCNParameters; ++i)
00595 if ((x[i] < fMCMCBoundaryMin[i]) || (x[i] > fMCMCBoundaryMax[i]))
00596 return false;
00597
00598 return true;
00599 }
00600
00601
00602 bool BCEngineMCMC::MCMCGetProposalPointMetropolis(int ichain, int ipar, std::vector<double> &x)
00603 {
00604
00605
00606 double proposal = MCMCTrialFunctionSingle(ichain, ipar);
00607
00608
00609 for (int i = 0; i < fMCMCNParameters; ++i)
00610 x[i] = fMCMCx[ichain * fMCMCNParameters + i];
00611
00612
00613 x[ipar] += proposal * (fMCMCBoundaryMax[ipar] - fMCMCBoundaryMin[ipar]);
00614
00615
00616 if ((x[ipar] < fMCMCBoundaryMin[ipar]) || (x[ipar] > fMCMCBoundaryMax[ipar]))
00617 return false;
00618
00619 return true;
00620 }
00621
00622
00623 bool BCEngineMCMC::MCMCGetNewPointMetropolis(int chain, int parameter)
00624 {
00625
00626 int index = chain * fMCMCNParameters;
00627
00628 fMCMCCurrentChain = chain;
00629
00630
00631 fMCMCNIterations[chain]++;
00632
00633
00634 if (!MCMCGetProposalPointMetropolis(chain, parameter, fMCMCxLocal))
00635 {
00636
00637 fMCMCNTrialsFalse[chain * fMCMCNParameters + parameter]++;
00638
00639
00640 MCMCCurrentPointInterface(fMCMCxLocal, chain, false);
00641
00642 return false;
00643 }
00644
00645
00646 double p0 = fMCMCprob[chain];
00647 double p1 = LogEval(fMCMCxLocal);
00648
00649
00650 bool accept = false;
00651
00652
00653 if (p1 >= p0)
00654 accept = true;
00655
00656 else
00657 {
00658 double r = log(fRandom->Rndm());
00659
00660 if(r < p1 - p0)
00661 accept = true;
00662 }
00663
00664
00665 if(accept)
00666 {
00667
00668 fMCMCNTrialsTrue[chain * fMCMCNParameters + parameter]++;
00669
00670
00671 for(int i = 0; i < fMCMCNParameters; ++i)
00672 {
00673
00674 fMCMCx[index + i] = fMCMCxLocal[i];
00675
00676
00677 fMCMCprob[chain] = p1;
00678 }
00679 }
00680 else
00681 {
00682
00683 fMCMCNTrialsFalse[chain * fMCMCNParameters + parameter]++;
00684 }
00685
00686
00687 MCMCCurrentPointInterface(fMCMCxLocal, chain, accept);
00688
00689 return accept;
00690 }
00691
00692
00693 bool BCEngineMCMC::MCMCGetNewPointMetropolis(int chain)
00694 {
00695
00696 int index = chain * fMCMCNParameters;
00697
00698 fMCMCCurrentChain = chain;
00699
00700
00701 fMCMCNIterations[chain]++;
00702
00703
00704 if (!MCMCGetProposalPointMetropolis(chain, fMCMCxLocal))
00705 {
00706
00707 for (int i = 0; i < fMCMCNParameters; ++i)
00708 fMCMCNTrialsFalse[chain * fMCMCNParameters + i]++;
00709
00710
00711 MCMCCurrentPointInterface(fMCMCxLocal, chain, false);
00712
00713 return false;
00714 }
00715
00716
00717 double p0 = fMCMCprob[chain];
00718 double p1 = LogEval(fMCMCxLocal);
00719
00720
00721 bool accept = false;
00722
00723
00724 if (p1 >= p0)
00725 accept = true;
00726
00727
00728 else
00729 {
00730 double r = log(fRandom->Rndm());
00731
00732 if(r < p1 - p0)
00733 accept = true;
00734 }
00735
00736
00737 if(accept)
00738 {
00739
00740 for (int i = 0; i < fMCMCNParameters; ++i)
00741 fMCMCNTrialsTrue[chain * fMCMCNParameters + i]++;
00742
00743
00744 for(int i = 0; i < fMCMCNParameters; ++i)
00745 {
00746
00747 fMCMCx[index + i] = fMCMCxLocal[i];
00748
00749
00750 fMCMCprob[chain] = p1;
00751 }
00752 }
00753 else
00754 {
00755
00756 for (int i = 0; i < fMCMCNParameters; ++i)
00757 fMCMCNTrialsFalse[chain * fMCMCNParameters + i]++;
00758 }
00759
00760
00761 MCMCCurrentPointInterface(fMCMCxLocal, chain, accept);
00762
00763 return accept;
00764 }
00765
00766
00767 void BCEngineMCMC::MCMCInChainCheckMaximum()
00768 {
00769
00770 for (int i = 0; i < fMCMCNChains; ++i)
00771 {
00772
00773 if (fMCMCprob[i] > fMCMCprobMax[i] || fMCMCNIterations[i] == 1)
00774 {
00775
00776 fMCMCprobMax[i] = fMCMCprob[i];
00777
00778
00779 for (int j = 0; j < fMCMCNParameters; ++j)
00780 fMCMCxMax[i * fMCMCNParameters + j] = fMCMCx[i * fMCMCNParameters + j];
00781 }
00782 }
00783 }
00784
00785
00786 void BCEngineMCMC::MCMCInChainUpdateStatistics()
00787 {
00788
00789 int npoints = fMCMCNTrialsTrue[0] + fMCMCNTrialsFalse[0];
00790
00791
00792 int nentries = fMCMCNParameters * fMCMCNChains;
00793
00794
00795 for (int i = 0; i < nentries; ++i) {
00796
00797 fMCMCxMean[i] += (fMCMCx[i] - fMCMCxMean[i]) / double(npoints);
00798
00799
00800 if (npoints > 1)
00801 fMCMCxVar[i] = (1.0 - 1./double(npoints)) * fMCMCxVar[i]
00802 + (fMCMCx[i] - fMCMCxMean[i]) * (fMCMCx[i] - fMCMCxMean[i]) / double(npoints - 1);
00803 }
00804
00805
00806 for (int i = 0; i < fMCMCNChains; ++i) {
00807
00808 fMCMCprobMean[i] += (fMCMCprob[i] - fMCMCprobMean[i]) / double(npoints);
00809
00810
00811 if (npoints > 1)
00812 fMCMCprobVar[i] = (1.0 - 1/double(npoints)) * fMCMCprobVar[i]
00813 + (fMCMCprob[i] - fMCMCprobMean[i]) * (fMCMCprob[i] - fMCMCprobMean[i]) / double(npoints - 1);
00814 }
00815
00816 }
00817
00818
00819 void BCEngineMCMC::MCMCInChainFillHistograms()
00820 {
00821
00822 if (!fMCMCFlagFillHistograms)
00823 return;
00824
00825
00826 for (int i = 0; i < fMCMCNChains; ++i)
00827 {
00828
00829 for (int j = 0; j < fMCMCNParameters; ++j)
00830 if (fMCMCFlagsFillHistograms.at(j))
00831 fMCMCH1Marginalized[j]->Fill(fMCMCx[i * fMCMCNParameters + j]);
00832
00833
00834 int counter = 0;
00835
00836 for (int j = 0; j < fMCMCNParameters; ++j)
00837 for (int k = 0; k < j; ++k)
00838 {
00839 if (fMCMCFlagsFillHistograms.at(j) && fMCMCFlagsFillHistograms.at(k))
00840 fMCMCH2Marginalized[counter]->Fill(fMCMCx[i*fMCMCNParameters+k],fMCMCx[i* fMCMCNParameters+j]);
00841 counter ++;
00842 }
00843 }
00844 }
00845
00846
00847 void BCEngineMCMC::MCMCInChainTestConvergenceAllChains()
00848 {
00849
00850 int npoints = fMCMCNTrialsTrue[0] + fMCMCNTrialsFalse[0];
00851
00852 if (fMCMCNChains > 1 && npoints > 1)
00853 {
00854
00855 bool flag_convergence = true;
00856
00857
00858 for (int iparameters = 0; iparameters < fMCMCNParameters; ++iparameters){
00859
00860 std::vector<double> means(fMCMCNChains);
00861 std::vector<double> variances(fMCMCNChains);
00862
00863
00864 for (int ichains = 0; ichains < fMCMCNChains; ++ichains) {
00865
00866 int index = ichains * fMCMCNParameters + iparameters;
00867 means[ichains] = fMCMCxMean[index];
00868 variances[ichains] = fMCMCxVar[index];
00869 }
00870
00871 fMCMCRValueParameters[iparameters] = BCMath::Rvalue(means, variances, npoints, fMCMCRValueUseStrict);
00872
00873
00874 if (! ((fMCMCRValueParameters[iparameters]-1.0) < fMCMCRValueParametersCriterion))
00875 flag_convergence = false;
00876
00877
00878 }
00879
00880 fMCMCRValue = BCMath::Rvalue(fMCMCprobMean, fMCMCprobVar, npoints, fMCMCRValueUseStrict);
00881
00882
00883 if (!((fMCMCRValue - 1.0) < fMCMCRValueCriterion))
00884 flag_convergence = false;
00885
00886
00887 if (fMCMCNIterationsConvergenceGlobal == -1 && flag_convergence == true)
00888 fMCMCNIterationsConvergenceGlobal = fMCMCNIterations[0] / fMCMCNParameters;
00889 }
00890 }
00891
00892 void BCEngineMCMC::MCMCInChainWriteChains()
00893 {
00894
00895 for (int i = 0; i < fMCMCNChains; ++i)
00896 fMCMCTrees[i]->Fill();
00897 }
00898
00899
00900 double BCEngineMCMC::LogEval(const std::vector<double> & )
00901 {
00902
00903
00904 return 0.0;
00905 }
00906
00907
00908 int BCEngineMCMC::MCMCMetropolisPreRun()
00909 {
00910
00911 BCLog::OutSummary("Pre-run Metropolis MCMC...");
00912
00913
00914 MCMCInitialize();
00915 MCMCInitializeMarkovChains();
00916
00917
00918 int ndigits = (int)log10(fMCMCNParameters) +1;
00919 if(ndigits<4)
00920 ndigits=4;
00921
00922
00923 MCMCResetRunStatistics();
00924 fMCMCNIterationsConvergenceGlobal = -1;
00925
00926
00927 BCLog::OutSummary(Form(" --> Perform MCMC pre-run with %i chains, each with maximum %i iterations", fMCMCNChains, fMCMCNIterationsMax));
00928
00929
00930
00931 bool tempflag_writetofile = fMCMCFlagWriteChainToFile;
00932 fMCMCFlagWriteChainToFile = false;
00933
00934
00935 fMCMCCurrentIteration = 1;
00936 int counterupdate = 1;
00937 bool convergence = false;
00938 bool flagefficiency = false;
00939
00940
00941 std::vector<double> efficiency;
00942 efficiency.assign(fMCMCNParameters * fMCMCNChains, 0.0);
00943
00944
00945
00946
00947
00948 int updateLimit = ( fMCMCNIterationsUpdateMax<fMCMCNIterationsUpdate*(fMCMCNParameters) && fMCMCNIterationsUpdateMax>0 ) ?
00949 fMCMCNIterationsUpdateMax : fMCMCNIterationsUpdate*(fMCMCNParameters);
00950
00951
00952 for (int ichains = 0; ichains < fMCMCNChains; ++ichains) {
00953
00954 for (int iparameter = 0; iparameter < fMCMCNParameters; ++iparameter){
00955
00956 int index = ichains * fMCMCNParameters + iparameter;
00957
00958 fMCMCNTrialsTrue[index] = 0;
00959 fMCMCNTrialsFalse[index] = 0;
00960 fMCMCxMean[index] = fMCMCx[index];
00961 fMCMCxVar[index] = 0;
00962 }
00963 fMCMCprobMean[ichains] = fMCMCprob[ichains];
00964 fMCMCprobVar[ichains] = 0;
00965 }
00966
00967
00968 fMCMCPhase = 1;
00969 fMCMCCycle = 0;
00970
00971
00972
00973
00974
00975
00976 while (fMCMCCurrentIteration < fMCMCNIterationsPreRunMin || (fMCMCCurrentIteration < fMCMCNIterationsMax && !(convergence && flagefficiency)))
00977 {
00978
00979
00980
00981
00982
00983 convergence = false;
00984
00985
00986 fMCMCNIterationsConvergenceGlobal = -1;
00987
00988
00989
00990
00991
00992
00993 if (fMCMCFlagOrderParameters)
00994 {
00995
00996 for (int iparameters = 0; iparameters < fMCMCNParameters; ++iparameters)
00997 {
00998
00999 for (int ichains = 0; ichains < fMCMCNChains; ++ichains)
01000 MCMCGetNewPointMetropolis(ichains, iparameters);
01001
01002
01003 MCMCInChainCheckMaximum();
01004 }
01005 }
01006
01007
01008 else
01009 {
01010
01011 for (int ichains = 0; ichains < fMCMCNChains; ++ichains)
01012
01013 MCMCGetNewPointMetropolis(ichains);
01014
01015
01016 MCMCInChainCheckMaximum();
01017 }
01018
01019
01020
01021
01022
01023
01024 if ( fMCMCCurrentIteration > 0 && fMCMCCurrentIteration % fMCMCNIterationsUpdate == 0 )
01025 BCLog::OutDetail(Form(" --> Iteration %i", fMCMCNIterations[0]/fMCMCNParameters));
01026
01027
01028
01029
01030
01031 if (counterupdate > 1)
01032 MCMCInChainUpdateStatistics();
01033
01034
01035
01036
01037
01038
01039
01040 if ( counterupdate % updateLimit == 0 && counterupdate > 0 && fMCMCCurrentIteration >= fMCMCNIterationsPreRunMin)
01041
01042 {
01043
01044
01045
01046
01047 bool rvalues_ok = true;
01048
01049 static bool has_converged = false;
01050
01051
01052
01053 fMCMCNIterationsConvergenceGlobal = -1;
01054
01055
01056
01057
01058
01059
01060 MCMCInChainTestConvergenceAllChains();
01061
01062
01063 if (fMCMCNIterationsConvergenceGlobal > 0)
01064 convergence = true;
01065
01066
01067 if (convergence && fMCMCNChains > 1)
01068 BCLog::OutDetail(Form(" * Convergence status: Set of %i Markov chains converged within %i iterations.", fMCMCNChains, fMCMCNIterationsConvergenceGlobal));
01069 else if (!convergence && fMCMCNChains > 1)
01070 {
01071 BCLog::OutDetail(Form(" * Convergence status: Set of %i Markov chains did not converge after %i iterations.", fMCMCNChains, fMCMCCurrentIteration));
01072
01073 BCLog::OutDetail(" - R-Values:");
01074 for (int iparameter = 0; iparameter < fMCMCNParameters; ++iparameter)
01075 {
01076 if(fabs(fMCMCRValueParameters[iparameter]-1.) < fMCMCRValueParametersCriterion)
01077 BCLog::OutDetail(Form( TString::Format(" parameter %%%di : %%.06f",ndigits), iparameter, fMCMCRValueParameters.at(iparameter)));
01078 else
01079 {
01080 BCLog::OutDetail(Form( TString::Format(" parameter %%%di : %%.06f <--",ndigits), iparameter, fMCMCRValueParameters.at(iparameter)));
01081 rvalues_ok = false;
01082 }
01083 }
01084 if(fabs(fMCMCRValue-1.) < fMCMCRValueCriterion)
01085 BCLog::OutDetail(Form(" log-likelihood : %.06f", fMCMCRValue));
01086 else
01087 {
01088 BCLog::OutDetail(Form(" log-likelihood : %.06f <--", fMCMCRValue));
01089 rvalues_ok = false;
01090 }
01091 }
01092
01093
01094 if(!has_converged)
01095 if(rvalues_ok)
01096 has_converged = true;
01097
01098
01099
01100
01101
01102
01103 flagefficiency = true;
01104
01105 bool flagprintefficiency = true;
01106
01107
01108 for (int ichains = 0; ichains < fMCMCNChains; ++ichains)
01109 {
01110
01111 for (int iparameter = 0; iparameter < fMCMCNParameters; ++iparameter)
01112 {
01113
01114 int index = ichains * fMCMCNParameters + iparameter;
01115
01116
01117 efficiency[index] = double(fMCMCNTrialsTrue[index]) / double(fMCMCNTrialsTrue[index] + fMCMCNTrialsFalse[index]);
01118
01119
01120 if (efficiency[index] < fMCMCEfficiencyMin && fMCMCTrialFunctionScaleFactor[index] > .01)
01121 {
01122 if (flagprintefficiency)
01123 {
01124 BCLog::OutDetail(Form(" * Efficiency status: Efficiencies not within pre-defined range."));
01125 BCLog::OutDetail(Form(" - Efficiencies:"));
01126 flagprintefficiency = false;
01127 }
01128
01129 double fscale=2.;
01130 if(has_converged && fMCMCEfficiencyMin/efficiency[index] > 2.)
01131 fscale = 4.;
01132 fMCMCTrialFunctionScaleFactor[index] /= fscale;
01133
01134 BCLog::OutDetail(Form(" Efficiency of parameter %i dropped below %.2lf%% (eps = %.2lf%%) in chain %i. Set scale to %.4g",
01135 iparameter, 100. * fMCMCEfficiencyMin, 100. * efficiency[index], ichains, fMCMCTrialFunctionScaleFactor[index]));
01136 }
01137
01138
01139 else if (efficiency[index] > fMCMCEfficiencyMax && fMCMCTrialFunctionScaleFactor[index] < 1.0)
01140 {
01141 if (flagprintefficiency)
01142 {
01143 BCLog::OutDetail(Form(" * Efficiency status: Efficiencies not within pre-defined ranges."));
01144 BCLog::OutDetail(Form(" - Efficiencies:"));
01145 flagprintefficiency = false;
01146 }
01147
01148 fMCMCTrialFunctionScaleFactor[index] *= 2.;
01149
01150 BCLog::OutDetail(Form(" Efficiency of parameter %i above %.2lf%% (eps = %.2lf%%) in chain %i. Set scale to %.4g",
01151 iparameter, 100.0 * fMCMCEfficiencyMax, 100.0 * efficiency[index], ichains, fMCMCTrialFunctionScaleFactor[index]));
01152 }
01153
01154
01155 if ((efficiency[index] < fMCMCEfficiencyMin && fMCMCTrialFunctionScaleFactor[index] > .01)
01156 || (efficiency[index] > fMCMCEfficiencyMax && fMCMCTrialFunctionScaleFactor[index] < 1.))
01157 flagefficiency = false;
01158 }
01159 }
01160
01161
01162 if (flagefficiency)
01163 BCLog::OutDetail(Form(" * Efficiency status: Efficiencies within pre-defined ranges."));
01164
01165
01166
01167
01168
01169 counterupdate = 0;
01170
01171
01172 for (int ichains = 0; ichains < fMCMCNChains; ++ichains) {
01173
01174 for (int iparameter = 0; iparameter < fMCMCNParameters; ++iparameter){
01175
01176 int index = ichains * fMCMCNParameters + iparameter;
01177
01178 fMCMCNTrialsTrue[index] = 0;
01179 fMCMCNTrialsFalse[index] = 0;
01180 fMCMCxMean[index] = fMCMCx[index];
01181 fMCMCxVar[index] = 0;
01182 }
01183 fMCMCprobMean[ichains] = fMCMCprob[ichains];
01184 fMCMCprobVar[ichains] = 0;
01185 }
01186 }
01187
01188
01189
01190
01191
01192
01193 if (fMCMCFlagWritePreRunToFile)
01194 MCMCInChainWriteChains();
01195
01196
01197
01198
01199
01200 if (counterupdate == 0 && fMCMCCurrentIteration > 1)
01201 fMCMCCycle++;
01202 fMCMCCurrentIteration++;
01203 counterupdate++;
01204
01205 }
01206
01207
01208
01209
01210
01211
01212 if (fMCMCCurrentIteration<updateLimit)
01213 {
01214 BCLog::OutWarning(" Convergence never checked !");
01215 BCLog::OutWarning(" Increase maximum number of iterations in the pre-run /MCMCSetNIterationsMax()/");
01216 BCLog::OutWarning(" or decrease maximum number of iterations for update /MCMCSetNIterationsUpdateMax()/");
01217 }
01218
01219
01220
01221
01222
01223
01224 if (fMCMCNIterationsConvergenceGlobal > 0)
01225 fMCMCFlagConvergenceGlobal = true;
01226 else
01227 fMCMCFlagConvergenceGlobal = false;
01228
01229
01230 if (fMCMCFlagConvergenceGlobal && fMCMCNChains > 1 && !flagefficiency)
01231 BCLog::OutSummary(Form(" --> Set of %i Markov chains converged within %i iterations but could not adjust scales.", fMCMCNChains, fMCMCNIterationsConvergenceGlobal));
01232
01233 else if (fMCMCFlagConvergenceGlobal && fMCMCNChains > 1 && flagefficiency)
01234 BCLog::OutSummary(Form(" --> Set of %i Markov chains converged within %i iterations and all scales are adjusted.", fMCMCNChains, fMCMCNIterationsConvergenceGlobal));
01235
01236 else if (!fMCMCFlagConvergenceGlobal && (fMCMCNChains > 1) && flagefficiency)
01237 BCLog::OutSummary(Form(" --> Set of %i Markov chains did not converge within %i iterations.", fMCMCNChains, fMCMCNIterationsMax));
01238
01239 else if (!fMCMCFlagConvergenceGlobal && (fMCMCNChains > 1) && !flagefficiency)
01240 BCLog::OutSummary(Form(" --> Set of %i Markov chains did not converge within %i iterations and could not adjust scales.", fMCMCNChains, fMCMCNIterationsMax));
01241
01242 else if(fMCMCNChains == 1)
01243 BCLog::OutSummary(" --> No convergence criterion for a single chain defined.");
01244
01245 else
01246 BCLog::OutSummary(" --> Only one Markov chain. No global convergence criterion defined.");
01247
01248 BCLog::OutSummary(Form(" --> Markov chains ran for %i iterations.", fMCMCCurrentIteration));
01249
01250
01251
01252 std::vector<double> efficiencies;
01253
01254 for (int i = 0; i < fMCMCNParameters; ++i)
01255 efficiencies.push_back(0.);
01256
01257 BCLog::OutDetail(" --> Average efficiencies:");
01258 for (int i = 0; i < fMCMCNParameters; ++i)
01259 {
01260 for (int j = 0; j < fMCMCNChains; ++j)
01261 efficiencies[i] += efficiency[j * fMCMCNParameters + i] / double(fMCMCNChains);
01262
01263 BCLog::OutDetail(Form(TString::Format(" --> parameter %%%di : %%.02f%%%%",ndigits), i, 100. * efficiencies[i]));
01264 }
01265
01266
01267
01268 std::vector<double> scalefactors;
01269
01270 for (int i = 0; i < fMCMCNParameters; ++i)
01271 scalefactors.push_back(0.0);
01272
01273 BCLog::OutDetail(" --> Average scale factors:");
01274 for (int i = 0; i < fMCMCNParameters; ++i)
01275 {
01276 for (int j = 0; j < fMCMCNChains; ++j)
01277 scalefactors[i] += fMCMCTrialFunctionScaleFactor[j * fMCMCNParameters + i] / double(fMCMCNChains);
01278
01279 BCLog::OutDetail(Form( TString::Format(" --> parameter %%%di : %%.02f%%%%",ndigits), i, 100. * scalefactors[i]));
01280 }
01281
01282
01283 fMCMCFlagWriteChainToFile = tempflag_writetofile;
01284
01285
01286 fMCMCCurrentIteration = -1;
01287
01288
01289 fMCMCCurrentChain = -1;
01290
01291
01292 return 1;
01293 }
01294
01295
01296 int BCEngineMCMC::MCMCMetropolis()
01297 {
01298
01299 if (fMCMCFlagPreRun)
01300 MCMCMetropolisPreRun();
01301 else
01302 BCLog::OutWarning("BCEngineMCMC::MCMCMetropolis. Not running prerun. This can cause trouble if the data have changed.");
01303
01304
01305 BCLog::OutSummary( "Run Metropolis MCMC...");
01306
01307
01308 MCMCResetRunStatistics();
01309
01310
01311 fMCMCPhase = 2;
01312 fMCMCCycle = 0;
01313
01314
01315 BCLog::OutSummary(Form(" --> Perform MCMC run with %i chains, each with %i iterations.", fMCMCNChains, fMCMCNIterationsRun));
01316
01317
01318
01319
01320
01321
01322
01323
01324
01325
01326
01327 int nwrite = fMCMCNIterationsRun/10;
01328 if(nwrite < 100)
01329 nwrite=100;
01330 else if(nwrite < 500)
01331 nwrite=1000;
01332 else if(nwrite < 10000)
01333 nwrite=1000;
01334 else
01335 nwrite=10000;
01336
01337
01338 for (fMCMCCurrentIteration = 1; fMCMCCurrentIteration <= fMCMCNIterationsRun; ++fMCMCCurrentIteration)
01339 {
01340 if ( (fMCMCCurrentIteration)%nwrite == 0 )
01341 BCLog::OutDetail(Form(" --> iteration number %i (%.2f%%)", fMCMCCurrentIteration, (double)(fMCMCCurrentIteration)/(double)fMCMCNIterationsRun*100.));
01342
01343
01344 if (fMCMCFlagOrderParameters)
01345 {
01346
01347 for (int iparameters = 0; iparameters < fMCMCNParameters; ++iparameters)
01348 {
01349
01350 for (int ichains = 0; ichains < fMCMCNChains; ++ichains)
01351 MCMCGetNewPointMetropolis(ichains, iparameters);
01352
01353
01354 fMCMCCurrentChain = -1;
01355
01356
01357 MCMCInChainCheckMaximum();
01358
01359 }
01360
01361
01362 if ( fMCMCCurrentIteration % fMCMCNLag == 0)
01363 {
01364
01365 MCMCInChainFillHistograms();
01366
01367
01368 if (fMCMCFlagWriteChainToFile)
01369 MCMCInChainWriteChains();
01370
01371
01372 MCMCIterationInterface();
01373 }
01374 }
01375
01376 else
01377 {
01378
01379 for (int ichains = 0; ichains < fMCMCNChains; ++ichains)
01380
01381 MCMCGetNewPointMetropolis(ichains);
01382
01383
01384 fMCMCCurrentChain = -1;
01385
01386
01387 MCMCInChainCheckMaximum();
01388
01389
01390 if (fMCMCCurrentIteration % fMCMCNLag == 0)
01391 {
01392
01393 MCMCInChainFillHistograms();
01394
01395
01396 if (fMCMCFlagWriteChainToFile)
01397 MCMCInChainWriteChains();
01398
01399
01400 MCMCIterationInterface();
01401 }
01402 }
01403
01404 }
01405
01406
01407 BCLog::OutSummary(Form(" --> Markov chains ran for %i iterations.", fMCMCNIterationsRun));
01408
01409
01410
01411
01412 double probmax = fMCMCprobMax.at(0);
01413 int probmaxindex = 0;
01414
01415
01416 for (int i = 1; i < fMCMCNChains; ++i) {
01417 if (fMCMCprobMax.at(i) > probmax)
01418 {
01419 probmax = fMCMCprobMax.at(i);
01420 probmaxindex = i;
01421 }
01422 }
01423
01424 BCLog::OutDetail(" --> Global mode from MCMC:");
01425 BCLog::OutDebug(Form(" --> Posterior value: %g", probmax));
01426 int ndigits = (int) log10(fMCMCNParameters);
01427 for (int i = 0; i < fMCMCNParameters; ++i)
01428 BCLog::OutDetail(Form( TString::Format(" --> parameter %%%di: %%.4g", ndigits+1),
01429 i, fMCMCxMax[probmaxindex * fMCMCNParameters + i]));
01430
01431
01432 fMCMCCurrentIteration = -1;
01433
01434
01435 fMCMCCurrentChain = -1;
01436
01437
01438 fMCMCFlagRun = true;
01439
01440 return 1;
01441 }
01442
01443
01444 void BCEngineMCMC::MCMCResetRunStatistics()
01445 {
01446 for (int j = 0; j < fMCMCNChains; ++j)
01447 {
01448 fMCMCNIterations[j] = 0;
01449 fMCMCNTrialsTrue[j] = 0;
01450 fMCMCNTrialsFalse[j] = 0;
01451 fMCMCprobMean[j] = 0;
01452 fMCMCprobVar[j] = 0;
01453
01454 for (int k = 0; k < fMCMCNParameters; ++k)
01455 {
01456 fMCMCNTrialsTrue[j * fMCMCNParameters + k] = 0;
01457 fMCMCNTrialsFalse[j * fMCMCNParameters + k] = 0;
01458 }
01459 }
01460
01461
01462 for (int i = 0; i < int(fMCMCH1Marginalized.size()); ++i)
01463 if (fMCMCH1Marginalized[i])
01464 fMCMCH1Marginalized[i]->Reset();
01465
01466 for (int i = 0; i < int(fMCMCH2Marginalized.size()); ++i)
01467 if (fMCMCH2Marginalized[i])
01468 fMCMCH2Marginalized[i]->Reset();
01469
01470 fMCMCRValue = 100;
01471 }
01472
01473
01474 int BCEngineMCMC::MCMCAddParameter(double min, double max)
01475 {
01476
01477 fMCMCBoundaryMin.push_back(min);
01478 fMCMCBoundaryMax.push_back(max);
01479
01480
01481 fMCMCFlagsFillHistograms.push_back(true);
01482
01483
01484 fMCMCNParameters++;
01485
01486
01487 return fMCMCNParameters;
01488 }
01489
01490
01491 void BCEngineMCMC::MCMCInitializeMarkovChains()
01492 {
01493
01494 std::vector<double> x0;
01495
01496 for (int j = 0; j < fMCMCNChains; ++j)
01497 {
01498 x0.clear();
01499 for (int i = 0; i < fMCMCNParameters; ++i)
01500 x0.push_back(fMCMCx[j * fMCMCNParameters + i]);
01501 fMCMCprob[j] = LogEval(x0);
01502 }
01503
01504 x0.clear();
01505 }
01506
01507
01508 int BCEngineMCMC::MCMCResetResults()
01509 {
01510
01511 fMCMCNIterations.clear();
01512 fMCMCNTrialsTrue.clear();
01513 fMCMCNTrialsFalse.clear();
01514 fMCMCTrialFunctionScaleFactor.clear();
01515 fMCMCprobMean.clear();
01516 fMCMCprobVar.clear();
01517 fMCMCxMean.clear();
01518 fMCMCxVar.clear();
01519 fMCMCx.clear();
01520 fMCMCprob.clear();
01521 fMCMCxMax.clear();
01522 fMCMCprobMax.clear();
01523 fMCMCNIterationsConvergenceGlobal = -1;
01524 fMCMCRValueParameters.clear();
01525
01526 for (int i = 0; i < int(fMCMCH1Marginalized.size()); ++i)
01527 if (fMCMCH1Marginalized[i])
01528 delete fMCMCH1Marginalized[i];
01529
01530 for (int i = 0; i < int(fMCMCH2Marginalized.size()); ++i)
01531 if (fMCMCH2Marginalized[i])
01532 delete fMCMCH2Marginalized[i];
01533
01534
01535 fMCMCH1Marginalized.clear();
01536 fMCMCH2Marginalized.clear();
01537
01538
01539 fMCMCFlagPreRun = true;
01540 fMCMCFlagRun = false;
01541 fMCMCFlagConvergenceGlobal = false;
01542
01543
01544 return 1;
01545 }
01546
01547
01548 int BCEngineMCMC::MCMCInitialize()
01549 {
01550
01551 MCMCResetResults();
01552
01553
01554 fMCMCNIterations.assign(fMCMCNChains, 0);
01555 fMCMCprobMean.assign(fMCMCNChains, 0);
01556 fMCMCprobVar.assign(fMCMCNChains, 0);
01557 fMCMCprob.assign(fMCMCNChains, -1.0);
01558 fMCMCprobMax.assign(fMCMCNChains, -1.0);
01559
01560 fMCMCNTrialsTrue.assign(fMCMCNChains * fMCMCNParameters, 0);
01561 fMCMCNTrialsFalse.assign(fMCMCNChains * fMCMCNParameters, 0);
01562 fMCMCxMax.assign(fMCMCNChains * fMCMCNParameters, 0.);
01563 fMCMCxMean.assign(fMCMCNChains * fMCMCNParameters, 0);
01564 fMCMCxVar.assign(fMCMCNChains * fMCMCNParameters, 0);
01565
01566 fMCMCRValueParameters.assign(fMCMCNParameters, 0.);
01567
01568 if (fMCMCTrialFunctionScaleFactorStart.size() == 0)
01569 fMCMCTrialFunctionScaleFactor.assign(fMCMCNChains * fMCMCNParameters, 1.0);
01570 else
01571 for (int i = 0; i < fMCMCNChains; ++i)
01572 for (int j = 0; j < fMCMCNParameters; ++j)
01573 fMCMCTrialFunctionScaleFactor.push_back(fMCMCTrialFunctionScaleFactorStart.at(j));
01574
01575
01576 if (fMCMCFlagInitialPosition == 2)
01577 {
01578
01579 bool flag = true;
01580
01581
01582 if (int(fMCMCInitialPosition.size()) != (fMCMCNChains * fMCMCNParameters))
01583 {
01584 BCLog::OutError("BCEngine::MCMCInitialize : Length of vector containing initial positions does not have required length.");
01585 flag = false;
01586 }
01587
01588
01589 if (flag)
01590 {
01591 for (int j = 0; j < fMCMCNChains; ++j)
01592 for (int i = 0; i < fMCMCNParameters; ++i)
01593 if (fMCMCInitialPosition[j * fMCMCNParameters + i] < fMCMCBoundaryMin[i] ||
01594 fMCMCInitialPosition[j * fMCMCNParameters + i] > fMCMCBoundaryMax[i])
01595 {
01596 BCLog::OutError("BCEngine::MCMCInitialize : Initial position out of boundaries.");
01597 flag = false;
01598 }
01599 }
01600
01601
01602 if (!flag)
01603 fMCMCFlagInitialPosition = 1;
01604 }
01605
01606 if (fMCMCFlagInitialPosition == 0)
01607 for (int j = 0; j < fMCMCNChains; ++j)
01608 for (int i = 0; i < fMCMCNParameters; ++i)
01609 fMCMCx.push_back(fMCMCBoundaryMin[i] + .5 * (fMCMCBoundaryMax[i] - fMCMCBoundaryMin[i]));
01610
01611 else if (fMCMCFlagInitialPosition == 2)
01612 {
01613 for (int j = 0; j < fMCMCNChains; ++j)
01614 for (int i = 0; i < fMCMCNParameters; ++i)
01615 fMCMCx.push_back(fMCMCInitialPosition.at(j * fMCMCNParameters + i));
01616 }
01617
01618 else
01619 for (int j = 0; j < fMCMCNChains; ++j)
01620 for (int i = 0; i < fMCMCNParameters; ++i)
01621 fMCMCx.push_back(fMCMCBoundaryMin[i] + fRandom->Rndm() * (fMCMCBoundaryMax[i] - fMCMCBoundaryMin[i]));
01622
01623
01624 fMCMCxLocal.clear();
01625 for (int i = 0; i < fMCMCNParameters; ++i)
01626 fMCMCxLocal.push_back(fMCMCx[i]);
01627
01628
01629 for(int i = 0; i < fMCMCNParameters; ++i)
01630 {
01631 TH1D * h1 = 0;
01632 if (fMCMCFlagsFillHistograms.at(i))
01633 h1 = new TH1D(TString::Format("h1_%d_parameter_%i", BCLog::GetHIndex() ,i), "",
01634 fMCMCH1NBins[i], fMCMCBoundaryMin[i], fMCMCBoundaryMax[i]);
01635 fMCMCH1Marginalized.push_back(h1);
01636 }
01637
01638 for(int i = 0; i < fMCMCNParameters; ++i)
01639 for (int k = 0; k < i; ++k)
01640 {
01641 TH2D * h2 = 0;
01642 if (fMCMCFlagsFillHistograms.at(i) && fMCMCFlagsFillHistograms.at(k))
01643 h2 = new TH2D(Form("h2_%d_parameters_%i_vs_%i", BCLog::GetHIndex(), i, k), "",
01644 fMCMCH1NBins[k], fMCMCBoundaryMin[k], fMCMCBoundaryMax[k],
01645 fMCMCH1NBins[i], fMCMCBoundaryMin[i], fMCMCBoundaryMax[i] );
01646 fMCMCH2Marginalized.push_back(h2);
01647 }
01648
01649 return 1;
01650 }
01651
01652
01653 int BCEngineMCMC::SetMarginalized(int index, TH1D * h)
01654 {
01655 if((int)fMCMCH1Marginalized.size()<=index)
01656 return 0;
01657
01658 if(h==0)
01659 return 0;
01660
01661 if((int)fMCMCH1Marginalized.size()==index)
01662 fMCMCH1Marginalized.push_back(h);
01663 else
01664 fMCMCH1Marginalized[index]=h;
01665
01666 return index;
01667 }
01668
01669
01670 int BCEngineMCMC::SetMarginalized(int index1, int index2, TH2D * h)
01671 {
01672 int counter = 0;
01673 int index = 0;
01674
01675
01676 for(int i = 0; i < fMCMCNParameters; i++)
01677 for (int j = 0; j < i; ++j)
01678 {
01679 if(j == index1 && i == index2)
01680 index = counter;
01681 counter++;
01682 }
01683
01684 if((int)fMCMCH2Marginalized.size()<=index)
01685 return 0;
01686
01687 if(h==0)
01688 return 0;
01689
01690 if((int)fMCMCH2Marginalized.size()==index)
01691 fMCMCH2Marginalized.push_back(h);
01692 else
01693 fMCMCH2Marginalized[index]=h;
01694
01695 return index;
01696 }
01697
01698
01699