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