00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "BCModel.h"
00011
00012 #include "BCDataPoint.h"
00013 #include "BCDataSet.h"
00014 #include "BCParameter.h"
00015 #include "BCH1D.h"
00016 #include "BCH2D.h"
00017 #include "BCGoFTest.h"
00018 #include "BCLog.h"
00019 #include "BCErrorCodes.h"
00020 #include "BCMath.h"
00021 #include "BCModelOutput.h"
00022
00023 #include <TROOT.h>
00024 #include <TNamed.h>
00025 #include <TCanvas.h>
00026 #include <TPostScript.h>
00027 #include <TDirectory.h>
00028 #include <TFile.h>
00029 #include <TKey.h>
00030 #include <TTree.h>
00031 #include <TMath.h>
00032 #include <TGraph.h>
00033 #include <TH1.h>
00034 #include <TH2D.h>
00035 #include <TH1I.h>
00036 #include <TRandom3.h>
00037 #include <TF1.h>
00038 #include <Math/QuantFuncMathCore.h>
00039
00040 #include <algorithm>
00041 #include <set>
00042
00043 #include <fstream>
00044 #include <iomanip>
00045
00046
00047 BCModel::BCModel(const char * name)
00048 : BCIntegrate()
00049 {
00050 fNormalization = -1.;
00051 fDataSet = 0;
00052 fParameterSet = new BCParameterSet;
00053
00054 fIndex = -1;
00055 fPValue = -1;
00056 fPValueNDoF = -1;
00057 fChi2NDoF = -1;
00058
00059 fName = (char *) name;
00060
00061 fDataPointUpperBoundaries = 0;
00062 fDataPointLowerBoundaries = 0;
00063
00064 fErrorBandXY = 0;
00065
00066 fGoFNChains = 5;
00067 fGoFNIterationsMax = 100000;
00068 fGoFNIterationsRun = 2000;
00069
00070 flag_discrete = false;
00071
00072 fPriorConstantAll = false;
00073 fPriorConstantValue = 0;
00074 }
00075
00076
00077 BCModel::BCModel()
00078 : BCIntegrate()
00079 {
00080 fNormalization = -1.;
00081 fDataSet = 0;
00082 fParameterSet = new BCParameterSet();
00083
00084 fIndex = -1;
00085 fPValue = -1;
00086 fPValueNDoF = -1;
00087 fChi2NDoF = -1;
00088
00089 fName = "model";
00090 fDataPointUpperBoundaries = 0;
00091 fDataPointLowerBoundaries = 0;
00092
00093 fGoFNChains = 5;
00094 fGoFNIterationsMax = 100000;
00095 fGoFNIterationsRun = 2000;
00096
00097 flag_discrete = false;
00098
00099 fPriorConstantAll = false;
00100 fPriorConstantValue = 0;
00101 }
00102
00103
00104 BCModel::BCModel(const BCModel & bcmodel) : BCIntegrate(bcmodel)
00105 {
00106 fIndex = bcmodel.fIndex;
00107 fName = bcmodel.fName;
00108 fModelAPriori = bcmodel.fModelAPriori;
00109 fModelAPosteriori = bcmodel.fModelAPosteriori;
00110 for (int i = 0; i < int(bcmodel.fParameterSet->size()); ++i) {
00111 if (bcmodel.fParameterSet->at(i)) {
00112 fParameterSet->push_back(new BCParameter(*(bcmodel.fParameterSet->at(i))));
00113
00114 }
00115 else
00116 fParameterSet->push_back(0);
00117 }
00118 if (fDataSet)
00119 fDataSet = bcmodel.fDataSet;
00120 else
00121 fDataSet = 0;
00122 if (bcmodel.fNDataPointsMinimum)
00123 fNDataPointsMinimum = bcmodel.fNDataPointsMinimum;
00124 else
00125 fNDataPointsMinimum = 0;
00126 if (bcmodel.fNDataPointsMaximum)
00127 fNDataPointsMaximum = bcmodel.fNDataPointsMaximum;
00128 else
00129 fNDataPointsMaximum = 0;
00130 fPValue = bcmodel.fPValue;
00131 fChi2NDoF = bcmodel.fChi2NDoF;
00132 fPValueNDoF = bcmodel.fPValueNDoF;
00133 flag_discrete = bcmodel.flag_discrete;
00134 fGoFNIterationsMax = bcmodel.fGoFNIterationsMax;
00135 fGoFNIterationsRun = bcmodel.fGoFNIterationsRun;
00136 fGoFNChains = bcmodel.fGoFNChains;
00137 for (int i = 0; i < int(bcmodel.fPriorContainer.size()); ++i) {
00138 if (bcmodel.fPriorContainer.at(i))
00139 fPriorContainer.push_back(new TNamed(*bcmodel.fPriorContainer.at(i)));
00140 else
00141 fPriorContainer.push_back(0);
00142 }
00143 fPriorConstantAll = bcmodel.fPriorConstantAll;
00144 fPriorConstantValue = bcmodel.fPriorConstantValue;
00145 fPriorContainerConstant = bcmodel.fPriorContainerConstant;
00146 fPriorContainerInterpolate = bcmodel.fPriorContainerInterpolate;
00147 fNormalization = bcmodel.fNormalization;
00148 }
00149
00150
00151 BCModel::~BCModel()
00152 {
00153 for (unsigned int i = 0; i < GetNParameters(); ++i)
00154 delete fPriorContainer[i];
00155 fPriorContainer.clear();
00156
00157 delete fParameterSet;
00158
00159 delete fDataPointLowerBoundaries;
00160 delete fDataPointUpperBoundaries;
00161 }
00162
00163
00164 BCModel & BCModel::operator = (const BCModel & bcmodel)
00165 {
00166 BCIntegrate::operator=(bcmodel);
00167 fIndex = bcmodel.fIndex;
00168 fName = bcmodel.fName;
00169 fModelAPriori = bcmodel.fModelAPriori;
00170 fModelAPosteriori = bcmodel.fModelAPosteriori;
00171 for (int i = 0; i < int(bcmodel.fParameterSet->size()); ++i) {
00172 if (bcmodel.fParameterSet->at(i)) {
00173 fParameterSet->push_back(new BCParameter(*(bcmodel.fParameterSet->at(i))));
00174
00175 }
00176 else
00177 fParameterSet->push_back(0);
00178 }
00179 if (fDataSet)
00180 fDataSet = bcmodel.fDataSet;
00181 else
00182 fDataSet = 0;
00183 if (bcmodel.fNDataPointsMinimum)
00184 fNDataPointsMinimum = bcmodel.fNDataPointsMinimum;
00185 else
00186 fNDataPointsMinimum = 0;
00187 if (bcmodel.fNDataPointsMaximum)
00188 fNDataPointsMaximum = bcmodel.fNDataPointsMaximum;
00189 else
00190 fNDataPointsMaximum = 0;
00191 fPValue = bcmodel.fPValue;
00192 fChi2NDoF = bcmodel.fChi2NDoF;
00193 fPValueNDoF = bcmodel.fPValueNDoF;
00194 flag_discrete = bcmodel.flag_discrete;
00195 fGoFNIterationsMax = bcmodel.fGoFNIterationsMax;
00196 fGoFNIterationsRun = bcmodel.fGoFNIterationsRun;
00197 fGoFNChains = bcmodel.fGoFNChains;
00198 for (int i = 0; i < int(bcmodel.fPriorContainer.size()); ++i) {
00199 if (bcmodel.fPriorContainer.at(i))
00200 fPriorContainer.push_back(new TNamed(*bcmodel.fPriorContainer.at(i)));
00201 else
00202 fPriorContainer.push_back(0);
00203 }
00204 fPriorConstantAll = bcmodel.fPriorConstantAll;
00205 fPriorConstantValue = bcmodel.fPriorConstantValue;
00206 fPriorContainerConstant = bcmodel.fPriorContainerConstant;
00207 fPriorContainerInterpolate = bcmodel.fPriorContainerInterpolate;
00208 fNormalization = bcmodel.fNormalization;
00209
00210
00211 return *this;
00212 }
00213
00214
00215 int BCModel::GetNDataPoints()
00216 {
00217 int npoints = 0;
00218 if (fDataSet)
00219 npoints = fDataSet->GetNDataPoints();
00220 else {
00221 BCLog::OutWarning("BCModel::GetNDataPoints() : No data set defined.");
00222 return ERROR_NOEVENTS;
00223 }
00224
00225 return npoints;
00226 }
00227
00228
00229 BCDataPoint * BCModel::GetDataPoint(unsigned int index)
00230 {
00231 if (fDataSet)
00232 return fDataSet->GetDataPoint(index);
00233
00234 BCLog::OutWarning("BCModel::GetDataPoint : No data set defined.");
00235 return 0;
00236 }
00237
00238
00239 BCParameter * BCModel::GetParameter(int index)
00240 {
00241 if (!fParameterSet)
00242 return 0;
00243
00244 if (index < 0 || index >= (int)GetNParameters()) {
00245 BCLog::OutWarning(
00246 Form("BCModel::GetParameter : Parameter index %d not within range.",index));
00247 return 0;
00248 }
00249
00250 return fParameterSet->at(index);
00251 }
00252
00253
00254 BCParameter * BCModel::GetParameter(const char * name)
00255 {
00256 if (!fParameterSet)
00257 return 0;
00258
00259 int index = -1;
00260 for (unsigned int i = 0; i < GetNParameters(); i++)
00261 if (name == GetParameter(i)->GetName())
00262 index = i;
00263
00264 if (index < 0) {
00265 BCLog::OutWarning(
00266 Form("BCModel::GetParameter : Model %s has no parameter named '%s'",
00267 GetName().data(), name));
00268 return 0;
00269 }
00270
00271 return GetParameter(index);
00272 }
00273
00274
00275 double BCModel::GetBestFitParameter(unsigned int index)
00276 {
00277 if(index >= GetNParameters()) {
00278 BCLog::OutError("BCModel::GetBestFitParameter : Parameter index out of range, returning -1e+111.");
00279 return -1e+111;
00280 }
00281
00282 if(fBestFitParameters.size()==0) {
00283 BCLog::OutError("BCModel::GetBestFitParameter : Mode finding not yet run, returning center of the range.");
00284 return (GetParameter(index)->GetUpperLimit() + GetParameter(index)->GetLowerLimit() ) / 2.;
00285 }
00286
00287 return fBestFitParameters[index];
00288 }
00289
00290
00291 double BCModel::GetBestFitParameterError(unsigned int index)
00292 {
00293 if(index >= GetNParameters()) {
00294 BCLog::OutError("BCModel::GetBestFitParameterError : Parameter index out of range, returning -1.");
00295 return -1;
00296 }
00297
00298 if(fBestFitParameterErrors.size()==0) {
00299 BCLog::OutError("BCModel::GetBestFitParameterError : Mode finding not yet run, returning -1.");
00300 return -1.;
00301 }
00302
00303 if(fBestFitParameterErrors[index]<0.)
00304 BCLog::OutWarning("BCModel::GetBestFitParameterError : Parameter error not available, returning -1.");
00305
00306 return fBestFitParameterErrors[index];
00307 }
00308
00309
00310 double BCModel::GetBestFitParameterMarginalized(unsigned int index)
00311 {
00312 if(index >= GetNParameters()) {
00313 BCLog::OutError("BCModel::GetBestFitParameterMarginalized : Parameter index out of range, returning -1e+111.");
00314 return -1e+111;
00315 }
00316
00317 if(fBestFitParametersMarginalized.size()==0) {
00318 BCLog::OutError("BCModel::GetBestFitParameterMarginalized : MCMC not yet run, returning center of the range.");
00319 return (GetParameter(index)->GetUpperLimit() + GetParameter(index)->GetLowerLimit() ) / 2.;
00320 }
00321
00322 return fBestFitParametersMarginalized[index];
00323 }
00324
00325
00326 void BCModel::SetNbins(const char * parname, int nbins)
00327 {
00328 BCParameter * p = GetParameter(parname);
00329 if (!p) {
00330 BCLog::OutWarning(
00331 Form("BCModel::SetNbins : Parameter '%s' not found so Nbins not set",parname));
00332 return;
00333 }
00334
00335 BCIntegrate::SetNbins(nbins, p->GetIndex());
00336 }
00337
00338
00339 std::vector<double> BCModel::GetErrorBand(double level)
00340 {
00341 std::vector<double> errorband;
00342
00343 if (!fErrorBandXY)
00344 return errorband;
00345
00346 int nx = fErrorBandXY->GetNbinsX();
00347 errorband.assign(nx, 0.);
00348
00349
00350 for (int ix = 1; ix <= nx; ix++) {
00351 TH1D * temphist = fErrorBandXY->ProjectionY("temphist", ix, ix);
00352
00353 int nprobSum = 1;
00354 double q[1];
00355 double probSum[1];
00356 probSum[0] = level;
00357
00358 temphist->GetQuantiles(nprobSum, q, probSum);
00359
00360 errorband[ix - 1] = q[0];
00361 }
00362
00363 return errorband;
00364 }
00365
00366
00367 TGraph * BCModel::GetErrorBandGraph(double level1, double level2)
00368 {
00369 if (!fErrorBandXY)
00370 return 0;
00371
00372
00373 int nx = fErrorBandXY->GetNbinsX();
00374
00375 TGraph * graph = new TGraph(2 * nx);
00376 graph->SetFillStyle(1001);
00377 graph->SetFillColor(kYellow);
00378
00379
00380 std::vector<double> ymin = GetErrorBand(level1);
00381 std::vector<double> ymax = GetErrorBand(level2);
00382
00383 for (int i = 0; i < nx; i++) {
00384 graph->SetPoint(i, fErrorBandXY->GetXaxis()->GetBinCenter(i + 1), ymin[i]);
00385 graph->SetPoint(nx + i, fErrorBandXY->GetXaxis()->GetBinCenter(nx - i), ymax[nx - i - 1]);
00386 }
00387
00388 return graph;
00389 }
00390
00391
00392 TH2D * BCModel::GetErrorBandXY_yellow(double level, int nsmooth)
00393 {
00394 if (!fErrorBandXY)
00395 return 0;
00396
00397 int nx = fErrorBandXY->GetNbinsX();
00398 int ny = fErrorBandXY->GetNbinsY();
00399
00400
00401 TH2D * hist_tempxy = (TH2D*) fErrorBandXY->Clone(
00402 TString::Format("%s_sub_%f.2", fErrorBandXY->GetName(), level));
00403 hist_tempxy->Reset();
00404 hist_tempxy->SetFillColor(kYellow);
00405
00406
00407 for (int ix = 1; ix < nx; ix++) {
00408 BCH1D * hist_temp = new BCH1D();
00409
00410 TH1D * hproj = fErrorBandXY->ProjectionY("temphist", ix, ix);
00411 if (nsmooth > 0)
00412 hproj->Smooth(nsmooth);
00413
00414 hist_temp->SetHistogram(hproj);
00415
00416 TH1D * hist_temp_yellow = hist_temp->GetSmallestIntervalHistogram(level);
00417
00418 for (int iy = 1; iy <= ny; ++iy)
00419 hist_tempxy->SetBinContent(ix, iy, hist_temp_yellow->GetBinContent(iy));
00420
00421 delete hist_temp_yellow;
00422 delete hist_temp;
00423 }
00424
00425 return hist_tempxy;
00426 }
00427
00428
00429 TGraph * BCModel::GetFitFunctionGraph(const std::vector<double> ¶meters)
00430 {
00431 if (!fErrorBandXY)
00432 return 0;
00433
00434
00435 int nx = fErrorBandXY->GetNbinsX();
00436 TGraph * graph = new TGraph(nx);
00437
00438
00439 for (int i = 0; i < nx; i++) {
00440 double x = fErrorBandXY->GetXaxis()->GetBinCenter(i + 1);
00441
00442 std::vector<double> xvec;
00443 xvec.push_back(x);
00444 double y = FitFunction(xvec, parameters);
00445 xvec.clear();
00446
00447 graph->SetPoint(i, x, y);
00448 }
00449
00450 return graph;
00451 }
00452
00453
00454 TGraph * BCModel::GetFitFunctionGraph(const std::vector<double> ¶meters, double xmin, double xmax, int n)
00455 {
00456
00457 TGraph * graph = new TGraph(n + 1);
00458
00459 double dx = (xmax - xmin) / (double) n;
00460
00461
00462 for (int i = 0; i <= n; i++) {
00463 double x = (double) i * dx + xmin;
00464 std::vector<double> xvec;
00465 xvec.push_back(x);
00466 double y = FitFunction(xvec, parameters);
00467
00468 xvec.clear();
00469
00470 graph->SetPoint(i, x, y);
00471 }
00472
00473 return graph;
00474 }
00475
00476
00477 bool BCModel::GetFlagBoundaries()
00478 {
00479 if (!fDataPointLowerBoundaries)
00480 return false;
00481
00482 if (!fDataPointUpperBoundaries)
00483 return false;
00484
00485 if (fDataPointLowerBoundaries->GetNValues() != fDataSet->GetDataPoint(0)->GetNValues())
00486 return false;
00487
00488 if (fDataPointUpperBoundaries->GetNValues() != fDataSet->GetDataPoint(0)->GetNValues())
00489 return false;
00490
00491 return true;
00492 }
00493
00494
00495 void BCModel::SetSingleDataPoint(BCDataPoint * datapoint)
00496 {
00497
00498 BCDataSet * dataset = new BCDataSet();
00499
00500
00501 dataset->AddDataPoint(datapoint);
00502
00503
00504 SetDataSet(dataset);
00505 }
00506
00507
00508 void BCModel::SetSingleDataPoint(BCDataSet * dataset, unsigned int index)
00509 {
00510 if (index > dataset->GetNDataPoints())
00511 return;
00512
00513 SetSingleDataPoint(dataset->GetDataPoint(index));
00514 }
00515
00516
00517 void BCModel::SetDataBoundaries(unsigned int index, double lowerboundary, double upperboundary, bool fixed)
00518 {
00519
00520 if (!fDataSet) {
00521 BCLog::OutError("BCModel::SetDataBoundaries : Need to define data set first.");
00522 return;
00523 }
00524
00525
00526 if (index > fDataSet->GetDataPoint(0)->GetNValues()) {
00527 BCLog::OutError("BCModel::SetDataBoundaries : Index out of range.");
00528 return;
00529 }
00530
00531
00532 if (!fDataPointLowerBoundaries)
00533 fDataPointLowerBoundaries = new BCDataPoint(fDataSet->GetDataPoint(0)->GetNValues());
00534
00535 if (!fDataPointUpperBoundaries)
00536 fDataPointUpperBoundaries = new BCDataPoint(fDataSet->GetDataPoint(0)->GetNValues());
00537
00538 if (fDataFixedValues.size() == 0)
00539 fDataFixedValues.assign(fDataSet->GetDataPoint(0)->GetNValues(), false);
00540
00541
00542 fDataPointLowerBoundaries->SetValue(index, lowerboundary);
00543 fDataPointUpperBoundaries->SetValue(index, upperboundary);
00544 fDataFixedValues[index] = fixed;
00545 }
00546
00547
00548 void BCModel::SetErrorBandContinuous(bool flag)
00549 {
00550 fErrorBandContinuous = flag;
00551
00552 if (flag)
00553 return;
00554
00555
00556 fErrorBandX.clear();
00557
00558
00559 for (unsigned int i = 0; i < fDataSet->GetNDataPoints(); ++i)
00560 fErrorBandX.push_back(fDataSet->GetDataPoint(i)->GetValue(fFitFunctionIndexX));
00561 }
00562
00563
00564 int BCModel::AddParameter(const char * name, double lowerlimit, double upperlimit)
00565 {
00566
00567 BCParameter * parameter = new BCParameter(name, lowerlimit, upperlimit);
00568
00569 int flag_ok = AddParameter(parameter);
00570 if (flag_ok)
00571 delete parameter;
00572
00573 return flag_ok;
00574 }
00575
00576
00577 int BCModel::AddParameter(BCParameter * parameter)
00578 {
00579
00580 if (!fParameterSet) {
00581 BCLog::OutError("BCModel::AddParameter : Parameter set does not exist");
00582 return ERROR_PARAMETERSETDOESNOTEXIST;
00583 }
00584
00585
00586 int flag_exists = 0;
00587 for (unsigned int i = 0; i < GetNParameters(); i++)
00588 if (CompareStrings(parameter->GetName().data(), GetParameter(i)->GetName().data()) == 0)
00589 flag_exists = -1;
00590
00591 if (flag_exists < 0) {
00592 BCLog::OutError(Form(
00593 "BCModel::AddParameter : Parameter with name %s exists already. ",
00594 parameter->GetName().data()));
00595 return ERROR_PARAMETEREXISTSALREADY;
00596 }
00597
00598
00599 unsigned int index = fParameterSet->size();
00600 parameter->SetIndex(index);
00601
00602
00603 fParameterSet->push_back(parameter);
00604
00605
00606 fPriorContainer.push_back(0);
00607
00608
00609 fPriorContainerInterpolate.push_back(false);
00610
00611
00612 fPriorContainerConstant.push_back(false);
00613
00614
00615 SetParameters(fParameterSet);
00616
00617
00618 ResetResults();
00619
00620 return 0;
00621 }
00622
00623
00624 double BCModel::LogProbabilityNN(const std::vector<double> ¶meters)
00625 {
00626
00627 double logprob = LogLikelihood(parameters);
00628
00629
00630 logprob += LogAPrioriProbability(parameters);
00631
00632 return logprob;
00633 }
00634
00635
00636 double BCModel::LogProbability(const std::vector<double> ¶meters)
00637 {
00638
00639 if (fNormalization < 0. || fNormalization == 0.) {
00640 BCLog::OutError("BCModel::LogProbability. Normalization not available or zero.");
00641 return 0.;
00642 }
00643
00644 return LogProbabilityNN(parameters) - log(fNormalization);
00645 }
00646
00647
00648 double BCModel::LogAPrioriProbability(const std::vector<double> ¶meters)
00649 {
00650 double logprob = 0;
00651
00652 if(!fPriorConstantAll) {
00653
00654
00655 int npar = GetNParameters();
00656
00657
00658 for (int i = 0; i < npar; ++i) {
00659 if (fPriorContainer[i]) {
00660
00661 TF1 * f = dynamic_cast<TF1*>(fPriorContainer[i]);
00662 TH1 * h = dynamic_cast<TH1*>(fPriorContainer[i]);
00663
00664 if (f)
00665 logprob += log(f->Eval(parameters[i]));
00666 else if (h) {
00667 if(fPriorContainerInterpolate[i])
00668 logprob += log(h->Interpolate(parameters[i]));
00669 else
00670 logprob += log(h->GetBinContent(h->FindBin(parameters[i])));
00671 }
00672 else
00673 BCLog::OutError(Form(
00674 "BCModel::LogAPrioriProbability : Prior for parameter %s "
00675 "is defined but not recodnized.",
00676 GetParameter(i)->GetName().c_str()));
00677 }
00678
00679 else if (!fPriorContainerConstant[i]) {
00680 BCLog::OutWarning(Form(
00681 "BCModel::LogAPrioriProbability: Prior for parameter %s "
00682 "is undefined. Using constant prior to proceed.",
00683 GetParameter(i)->GetName().c_str()));
00684 logprob -= log(GetParameter(i)->GetRangeWidth());
00685 }
00686 }
00687 }
00688
00689
00690 logprob += fPriorConstantValue;
00691
00692 return logprob;
00693 }
00694
00695
00696 double BCModel::LogEval(const std::vector<double> ¶meters)
00697 {
00698 return LogProbabilityNN(parameters);
00699 }
00700
00701
00702 double BCModel::EvalSampling(const std::vector<double> ¶meters)
00703 {
00704 return SamplingFunction(parameters);
00705 }
00706
00707
00708 double BCModel::SamplingFunction(const std::vector<double> & )
00709 {
00710 double probability = 1.;
00711 for (std::vector<BCParameter *>::const_iterator it = fParameterSet->begin(); it != fParameterSet->end(); ++it)
00712 probability *= 1. / ((*it)->GetUpperLimit() - (*it)->GetLowerLimit());
00713 return probability;
00714 }
00715
00716
00717 double BCModel::Normalize()
00718 {
00719 if(fParameterSet->size()<1) {
00720 BCLog::OutError(Form("Normalize : No parameters defined in model \'%s\'. Aborting.",GetName().data()));
00721 return -1.;
00722 }
00723
00724 BCLog::OutSummary(Form("Model \'%s\': Normalizing probability",GetName().data()));
00725
00726 unsigned int n = GetNvar();
00727
00728
00729 if (n == 0) {
00730 SetParameters(fParameterSet);
00731 n = GetNvar();
00732 }
00733
00734
00735
00736 fNormalization = Integrate();
00737
00738 BCLog::OutDetail(Form(" --> Normalization factor : %.6g", fNormalization));
00739
00740 return fNormalization;
00741 }
00742
00743
00744 int BCModel::CheckParameters(const std::vector<double> ¶meters)
00745 {
00746
00747 if (!parameters.size() == fParameterSet->size())
00748 return ERROR_INVALIDNUMBEROFPARAMETERS;
00749
00750
00751 for (unsigned int i = 0; i < fParameterSet->size(); i++) {
00752 BCParameter * modelparameter = fParameterSet->at(i);
00753
00754 if (modelparameter->GetLowerLimit() > parameters.at(i)
00755 || modelparameter->GetUpperLimit() < parameters.at(i)) {
00756 BCLog::OutError(Form(
00757 "BCModel::CheckParameters : Parameter %s not within limits.",
00758 fParameterSet-> at(i)-> GetName().data()));
00759 return ERROR_PARAMETERNOTWITHINRANGE;
00760 }
00761 }
00762
00763 return 0;
00764 }
00765
00766
00767 void BCModel::FindMode(std::vector<double> start)
00768 {
00769 if(fParameterSet->size()<1) {
00770 BCLog::OutError(Form("FindMode : No parameters defined in model \'%s\'. Aborting.",GetName().data()));
00771 return;
00772 }
00773
00774
00775
00776 BCLog::OutSummary(Form("Model \'%s\': Finding mode using %s",GetName().data(), DumpOptimizationMethod().c_str()));
00777
00778
00779 SetParameters(fParameterSet);
00780
00781 switch (GetOptimizationMethod()) {
00782 case BCIntegrate::kOptSA:
00783 FindModeSA(start);
00784 return;
00785
00786 case BCIntegrate::kOptMinuit: {
00787 int printlevel = -1;
00788 if (BCLog::GetLogLevelScreen() <= BCLog::detail)
00789 printlevel = 0;
00790 if (BCLog::GetLogLevelScreen() <= BCLog::debug)
00791 printlevel = 1;
00792 BCIntegrate::FindModeMinuit(start, printlevel);
00793 return;
00794 }
00795
00796 case BCIntegrate::kOptMetropolis:
00797 MarginalizeAll();
00798 return;
00799
00800 default:
00801 BCLog::OutError(Form("BCModel::FindMode : Invalid mode finding method: %d",GetOptimizationMethod()));
00802 break;
00803 }
00804
00805 return;
00806 }
00807
00808
00809 void BCModel::FindModeMinuit(std::vector<double> start, int printlevel)
00810 {
00811 if(fParameterSet->size()<1) {
00812 BCLog::OutError(Form("FindModeMinuit : No parameters defined in model \'%s\'. Aborting.",GetName().data()));
00813 return;
00814 }
00815
00816
00817 SetParameters(fParameterSet);
00818 BCIntegrate::FindModeMinuit(start, printlevel);
00819 }
00820
00821
00822 void BCModel::WriteMode(const char * file)
00823 {
00824 ofstream ofi(file);
00825 if (!ofi.is_open()) {
00826 std::cerr << "Couldn't open file " << file << std::endl;
00827 return;
00828 }
00829
00830 int npar = fParameterSet->size();
00831 for (int i = 0; i < npar; i++)
00832 ofi << fBestFitParameters.at(i) << std::endl;
00833
00834 ofi << std::endl;
00835 ofi << "#######################################################################"
00836 << std::endl;
00837 ofi << "#" << std::endl;
00838 ofi << "# This file was created automatically by BCModel::WriteMode() call."
00839 << std::endl;
00840 ofi << "# It can be read in by call to BCModel::ReadMode()." << std::endl;
00841 ofi << "# Do not modify it unless you know what you're doing." << std::endl;
00842 ofi << "#" << std::endl;
00843 ofi << "#######################################################################"
00844 << std::endl;
00845 ofi << "#" << std::endl;
00846 ofi << "# Best fit parameters (mode) for model:" << std::endl;
00847 ofi << "# \'" << fName.data() << "\'" << std::endl;
00848 ofi << "#" << std::endl;
00849 ofi << "# Number of parameters: " << npar << std::endl;
00850 ofi << "# Parameters ordered as above:" << std::endl;
00851
00852 for (int i = 0; i < npar; i++) {
00853 ofi << "# " << i << ": ";
00854 ofi << fParameterSet->at(i)->GetName().data() << " = ";
00855 ofi << fBestFitParameters.at(i) << std::endl;
00856 }
00857
00858 ofi << "#" << std::endl;
00859 ofi << "########################################################################"
00860 << std::endl;
00861 }
00862
00863
00864 int BCModel::ReadMode(const char * file)
00865 {
00866 ifstream ifi(file);
00867 if (!ifi.is_open()) {
00868 BCLog::OutError(Form("BCModel::ReadMode : Couldn't open file %s.", file));
00869 return 0;
00870 }
00871
00872 int npar = fParameterSet->size();
00873 std::vector<double> mode;
00874
00875 int i = 0;
00876 while (i < npar && !ifi.eof()) {
00877 double a;
00878 ifi >> a;
00879 mode.push_back(a);
00880 i++;
00881 }
00882
00883 if (i < npar) {
00884 BCLog::OutError(Form("BCModel::ReadMode : Couldn't read mode from file %s.", file));
00885 BCLog::OutError(Form("BCModel::ReadMode : Expected %d parameters, found %d.", npar, i));
00886 return 0;
00887 }
00888
00889 BCLog::OutSummary(
00890 Form("# Read in best fit parameters (mode) for model \'%s\' from file %s:",
00891 fName.data(), file));
00892 SetMode(mode);
00893 for (int j = 0; j < npar; j++)
00894 BCLog::OutSummary(Form("# ->Parameter %d : %s = %e", j,
00895 fParameterSet->at(j)->GetName().data(), fBestFitParameters[j]));
00896
00897 BCLog::OutWarning("# ! Best fit values obtained before this call will be overwritten !");
00898
00899 return npar;
00900 }
00901
00902
00903 int BCModel::MarginalizeAll()
00904 {
00905 if(fParameterSet->size()<1) {
00906 BCLog::OutError(Form("MarginalizeAll : No parameters defined in model \'%s\'. Aborting.",GetName().data()));
00907 return 0;
00908 }
00909
00910 BCLog::OutSummary(Form("Running MCMC for model \'%s\'",GetName().data()));
00911
00912
00913 double dx = 0.;
00914 double dy = 0.;
00915
00916 if (fFitFunctionIndexX >= 0) {
00917 dx = (fDataPointUpperBoundaries->GetValue(fFitFunctionIndexX)
00918 - fDataPointLowerBoundaries->GetValue(fFitFunctionIndexX))
00919 / (double) fErrorBandNbinsX;
00920
00921 dy = (fDataPointUpperBoundaries->GetValue(fFitFunctionIndexY)
00922 - fDataPointLowerBoundaries->GetValue(fFitFunctionIndexY))
00923 / (double) fErrorBandNbinsY;
00924
00925 fErrorBandXY
00926 = new TH2D(TString::Format("errorbandxy_%d", BCLog::GetHIndex()), "",
00927 fErrorBandNbinsX,
00928 fDataPointLowerBoundaries->GetValue(fFitFunctionIndexX) - .5 * dx,
00929 fDataPointUpperBoundaries->GetValue(fFitFunctionIndexX) + .5 * dx,
00930 fErrorBandNbinsY,
00931 fDataPointLowerBoundaries->GetValue(fFitFunctionIndexY) - .5 * dy,
00932 fDataPointUpperBoundaries->GetValue(fFitFunctionIndexY) + .5 * dy);
00933 fErrorBandXY->SetStats(kFALSE);
00934
00935 for (int ix = 1; ix <= fErrorBandNbinsX; ++ix)
00936 for (int iy = 1; iy <= fErrorBandNbinsX; ++iy)
00937 fErrorBandXY->SetBinContent(ix, iy, 0.);
00938 }
00939
00940
00941 MCMCMetropolis();
00942
00943
00944 StoreMode();
00945
00946 return 1;
00947 }
00948
00949
00950 BCH1D * BCModel::GetMarginalized(BCParameter * parameter)
00951 {
00952 if (!parameter) {
00953
00954
00955 return 0;
00956 }
00957
00958 int index = parameter->GetIndex();
00959 if (fMCMCFlagsFillHistograms[index] == false) {
00960
00961
00962
00963
00964 return 0;
00965 }
00966
00967
00968 TH1D * hist = MCMCGetH1Marginalized(index);
00969 if (!hist)
00970 return 0;
00971
00972
00973 hist->SetName(Form("hist_%s_%s", GetName().data(), parameter->GetName().data()));
00974 hist->SetXTitle(parameter->GetName().data());
00975 hist->SetYTitle(Form("p(%s|data)", parameter->GetName().data()));
00976 hist->SetStats(kFALSE);
00977
00978
00979 BCH1D * hprob = new BCH1D();
00980 hprob->SetHistogram(hist);
00981
00982
00983 double bestfit = hprob->GetMode();
00984
00985 if (fBestFitParametersMarginalized.empty())
00986 fBestFitParametersMarginalized.assign(GetNParameters(), 0.0);
00987
00988 fBestFitParametersMarginalized[index] = bestfit;
00989
00990 hprob->SetGlobalMode(fBestFitParameters.at(index));
00991
00992 return hprob;
00993 }
00994
00995
00996 int BCModel::ReadMarginalizedFromFile(const char * file)
00997 {
00998 TFile * froot = new TFile(file);
00999 if (!froot->IsOpen()) {
01000 BCLog::OutError(Form("BCModel::ReadMarginalizedFromFile : Couldn't open file %s.", file));
01001 return 0;
01002 }
01003
01004
01005
01006
01007
01008
01009
01010 MCMCInitialize();
01011
01012 int k = 0;
01013 int n = GetNParameters();
01014 for (int i = 0; i < n; i++) {
01015 BCParameter * a = GetParameter(i);
01016 TKey * key = froot->GetKey(TString::Format("hist_%s_%s",GetName().data(), a->GetName().data()));
01017 if (key) {
01018 TH1D * h1 = (TH1D*) key->ReadObjectAny(TH1D::Class());
01019 h1->SetDirectory(0);
01020 if (SetMarginalized(i, h1))
01021 k++;
01022 }
01023 else
01024 BCLog::OutWarning(
01025 Form("BCModel::ReadMarginalizedFromFile : Couldn't read histogram \"hist_%s_%s\" from file %s.",
01026 GetName().data(), a->GetName().data(), file));
01027 }
01028
01029 for (int i = 0; i < n - 1; i++) {
01030 for (int j = i + 1; j < n; j++) {
01031 BCParameter * a = GetParameter(i);
01032 BCParameter * b = GetParameter(j);
01033 TKey * key = froot->GetKey(
01034 TString::Format("hist_%s_%s_%s",GetName().data(), a->GetName().data(),b->GetName().data()));
01035 if (key) {
01036 TH2D * h2 = (TH2D*) key->ReadObjectAny(TH2D::Class());
01037 h2->SetDirectory(0);
01038 if (SetMarginalized(i, j, h2))
01039 k++;
01040 }
01041 else
01042 BCLog::OutWarning(
01043 Form("BCModel::ReadMarginalizedFromFile : Couldn't read histogram \"hist_%s_%s_%s\" from file %s.",
01044 GetName().data(), a->GetName().data(), b->GetName().data(), file));
01045 }
01046 }
01047
01048 froot->Close();
01049
01050 return k;
01051 }
01052
01053
01054 int BCModel::ReadErrorBandFromFile(const char * file)
01055 {
01056 TFile * froot = new TFile(file);
01057 if (!froot->IsOpen()) {
01058 BCLog::OutError(Form("BCModel::ReadErrorBandFromFile. Couldn't open file %s.", file));
01059 return 0;
01060 }
01061
01062 int r = 0;
01063
01064 TH2D * h2 = (TH2D*) froot->Get("errorbandxy");
01065 if (h2) {
01066 h2->SetDirectory(0);
01067 h2->SetName(TString::Format("errorbandxy_%d", BCLog::GetHIndex()));
01068 SetErrorBandHisto(h2);
01069 r = 1;
01070 }
01071 else
01072 BCLog::OutWarning(
01073 Form("BCModel::ReadErrorBandFromFile : Couldn't read histogram \"errorbandxy\" from file %s.",file));
01074
01075 froot->Close();
01076
01077 return r;
01078 }
01079
01080
01081 int BCModel::PrintAllMarginalized1D(const char * filebase)
01082 {
01083 if (fMCMCH1Marginalized.size() == 0) {
01084 BCLog::OutError("BCModel::PrintAllMarginalized : Marginalized distributions not available.");
01085 return 0;
01086 }
01087
01088 int n = GetNParameters();
01089 for (int i = 0; i < n; i++) {
01090 BCParameter * a = GetParameter(i);
01091 if (GetMarginalized(a))
01092 GetMarginalized(a)->Print(Form("%s_1D_%s.ps", filebase, a->GetName().data()));
01093 }
01094
01095 return n;
01096 }
01097
01098
01099 int BCModel::PrintAllMarginalized2D(const char * filebase)
01100 {
01101 if (fMCMCH2Marginalized.size() == 0) {
01102 BCLog::OutError("BCModel::PrintAllMarginalized : Marginalized distributions not available.");
01103 return 0;
01104 }
01105
01106 int k = 0;
01107 int n = GetNParameters();
01108 for (int i = 0; i < n - 1; i++) {
01109 for (int j = i + 1; j < n; j++) {
01110 BCParameter * a = GetParameter(i);
01111 BCParameter * b = GetParameter(j);
01112
01113 double meana = (a->GetLowerLimit() + a->GetUpperLimit()) / 2.;
01114 double deltaa = (a->GetUpperLimit() - a->GetLowerLimit());
01115 if (deltaa <= 1e-7 * meana)
01116 continue;
01117
01118 double meanb = (b->GetLowerLimit() + b->GetUpperLimit()) / 2.;
01119 double deltab = (b->GetUpperLimit() - b->GetLowerLimit());
01120 if (deltab <= 1e-7 * meanb)
01121 continue;
01122
01123 if (GetMarginalized(a, b))
01124 GetMarginalized(a, b)->Print(
01125 Form("%s_2D_%s_%s.ps",filebase, a->GetName().data(), b->GetName().data()));
01126 k++;
01127 }
01128 }
01129
01130 return k;
01131 }
01132
01133
01134 int BCModel::PrintAllMarginalized(const char * file, unsigned int hdiv, unsigned int vdiv)
01135 {
01136 if (!fMCMCFlagFillHistograms) {
01137 BCLog::OutError("BCModel::PrintAllMarginalized : Marginalized distributions not filled.");
01138 return 0;
01139 }
01140
01141 int npar = GetNParameters();
01142
01143 if (fMCMCH1Marginalized.size() == 0 || (fMCMCH2Marginalized.size() == 0
01144 && npar > 1)) {
01145 BCLog::OutError(
01146 "BCModel::PrintAllMarginalized : Marginalized distributions not available.");
01147 return 0;
01148 }
01149
01150
01151 if (fMCMCH1Marginalized.size() == 1 && fMCMCH2Marginalized.size() == 0) {
01152 BCParameter * a = GetParameter(0);
01153 if (GetMarginalized(a))
01154 GetMarginalized(a)->Print(file);
01155 return 1;
01156 }
01157
01158 int c_width = 600;
01159 int c_height = 600;
01160
01161 int type = 112;
01162
01163 if (hdiv > vdiv) {
01164 if (hdiv > 3) {
01165 c_width = 1000;
01166 c_height = 700;
01167 }
01168 else {
01169 c_width = 800;
01170 c_height = 600;
01171 }
01172 }
01173 else if (hdiv < vdiv) {
01174 if (hdiv > 3) {
01175 c_height = 1000;
01176 c_width = 700;
01177 }
01178 else {
01179 c_height = 800;
01180 c_width = 600;
01181 }
01182 type = 111;
01183 }
01184
01185
01186 int nplots2d = npar * (npar - 1) / 2;
01187 int nplots = npar + nplots2d;
01188
01189
01190 BCLog::OutSummary(
01191 Form(
01192 "Printing all marginalized distributions (%d x 1D + %d x 2D = %d) into file %s",
01193 npar, nplots2d, nplots, file));
01194 if (nplots > 100)
01195 BCLog::OutDetail("This can take a while...");
01196
01197
01198 TCanvas * c = new TCanvas("c", "canvas", c_width, c_height);
01199
01200 TPostScript * ps = new TPostScript(file, type);
01201
01202 if (type == 112)
01203 ps->Range(24, 16);
01204 else
01205 ps->Range(16, 24);
01206
01207
01208 ps->NewPage();
01209 c->cd();
01210 c->Clear();
01211 c->Divide(hdiv, vdiv);
01212
01213 int n = 0;
01214 for (int i = 0; i < npar; i++) {
01215
01216 BCParameter * a = GetParameter(i);
01217
01218
01219 if (!GetMarginalized(a))
01220 continue;
01221
01222
01223 if (GetMarginalized(a)->GetHistogram()->Integral() <= 0)
01224 continue;
01225
01226
01227 if (i != 0 && i % (hdiv * vdiv) == 0) {
01228 c->Update();
01229 ps->NewPage();
01230 c->cd();
01231 c->Clear();
01232 c->Divide(hdiv, vdiv);
01233 }
01234
01235
01236 c->cd(i % (hdiv * vdiv) + 1);
01237
01238
01239 if (a->GetRangeWidth() == 0)
01240 GetMarginalized(a)->Draw(4, a->GetLowerLimit());
01241 else
01242 GetMarginalized(a)->Draw();
01243
01244 n++;
01245
01246 if (n % 100 == 0)
01247 BCLog::OutDetail(Form(" --> %d plots done", n));
01248 }
01249
01250 if (n > 0) {
01251 c->Update();
01252 }
01253
01254
01255 int k = 0;
01256 for (int i = 0; i < npar - 1; i++) {
01257 for (int j = i + 1; j < npar; j++) {
01258
01259 BCParameter * a = GetParameter(i);
01260 BCParameter * b = GetParameter(j);
01261
01262
01263 if (!GetMarginalized(a, b))
01264 continue;
01265
01266
01267 if (GetMarginalized(a, b)->GetHistogram()->Integral() <= 0)
01268 continue;
01269
01270
01271 if ((k != 0 && k % (hdiv * vdiv) == 0) || k == 0) {
01272 c->Update();
01273 ps->NewPage();
01274 c->cd();
01275 c->Clear();
01276 c->Divide(hdiv, vdiv);
01277 }
01278
01279
01280 c->cd(k % (hdiv * vdiv) + 1);
01281
01282 double meana = (a->GetLowerLimit() + a->GetUpperLimit()) / 2.;
01283 double deltaa = (a->GetUpperLimit() - a->GetLowerLimit());
01284 if (deltaa <= 1e-7 * meana)
01285 continue;
01286
01287 double meanb = (b->GetLowerLimit() + b->GetUpperLimit()) / 2.;
01288 double deltab = (b->GetUpperLimit() - b->GetLowerLimit());
01289 if (deltab <= 1e-7 * meanb)
01290 continue;
01291
01292 GetMarginalized(a, b)->Draw(52);
01293 k++;
01294
01295 if ((n + k) % 100 == 0)
01296 BCLog::OutDetail(Form(" --> %d plots done", n + k));
01297 }
01298 }
01299
01300 if ((n + k) > 100 && (n + k) % 100 != 0)
01301 BCLog::OutDetail(Form(" --> %d plots done", n + k));
01302
01303 c->Update();
01304 ps->Close();
01305
01306 delete c;
01307 delete ps;
01308
01309
01310 return n + k;
01311 }
01312
01313
01314 BCH2D * BCModel::GetMarginalized(BCParameter * par1, BCParameter * par2)
01315 {
01316 int index1 = par1->GetIndex();
01317 int index2 = par2->GetIndex();
01318
01319 if (fMCMCFlagsFillHistograms[index1] == false || fMCMCFlagsFillHistograms[index2] == false) {
01320
01321
01322
01323
01324 return 0;
01325 }
01326
01327 if (index1 == index2) {
01328
01329
01330 return 0;
01331 }
01332
01333 if (par1->GetRangeWidth() == 0 || par2->GetRangeWidth() == 0){
01334 return 0;
01335 }
01336
01337 BCParameter * npar1 = par1;
01338 BCParameter * npar2 = par2;
01339
01340 if (index1 > index2) {
01341 npar1 = par2;
01342 npar2 = par1;
01343
01344 int itmp = index1;
01345 index1 = index2;
01346 index2 = itmp;
01347 }
01348
01349
01350 TH2D * hist = MCMCGetH2Marginalized(index1, index2);
01351
01352 if (hist == 0)
01353 return 0;
01354
01355 BCH2D * hprob = new BCH2D();
01356
01357
01358 hist->SetName(Form("hist_%s_%s_%s", GetName().data(), npar1->GetName().data(), npar2->GetName().data()));
01359 hist->SetXTitle(Form("%s", npar1->GetName().data()));
01360 hist->SetYTitle(Form("%s", npar2->GetName().data()));
01361 hist->SetStats(kFALSE);
01362
01363 double gmode[] = { fBestFitParameters.at(index1), fBestFitParameters.at(index2) };
01364 hprob->SetGlobalMode(gmode);
01365
01366
01367 hprob->SetHistogram(hist);
01368
01369 return hprob;
01370 }
01371
01372
01373 double BCModel::GetPvalueFromChi2(const std::vector<double> &par, int sigma_index)
01374 {
01375 double ll = LogLikelihood(par);
01376 int n = GetNDataPoints();
01377
01378 double sum_sigma = 0;
01379 for (int i = 0; i < n; i++)
01380 sum_sigma += log(GetDataPoint(i)->GetValue(sigma_index));
01381
01382 double chi2 = -2. * (ll + (double) n / 2. * log(2. * M_PI) + sum_sigma);
01383
01384 fPValue = TMath::Prob(chi2, n);
01385
01386 return fPValue;
01387 }
01388
01389
01390 std::vector<double> BCModel::GetChi2Runs(int , int )
01391 {
01392 std::vector<double> x;
01393 return x;
01394 }
01395
01396
01397 double BCModel::GetPvalueFromChi2Johnson(std::vector<double> par)
01398 {
01399 double chi2 = GetChi2Johnson(par);
01400
01401 fPValue = TMath::Prob(chi2, NumberBins() - 1);
01402 return fPValue;
01403 }
01404
01405
01406 double BCModel::GetChi2Johnson(std::vector<double> par, int nBins)
01407 {
01408 typedef unsigned int uint;
01409
01410
01411 int n = GetNDataPoints();
01412
01413 if (nBins < 0)
01414 nBins = NumberBins();
01415
01416
01417 std::vector<double> a;
01418 for (int i = 0; i <= nBins; i++)
01419 a.push_back(i / double(nBins));
01420
01421
01422 TH1I * hist = new TH1I("h1", "h1 title", nBins, 0., 1.);
01423
01424
01425
01426 if (flag_discrete) {
01427
01428 for (int j = 0; j < n; j++) {
01429
01430 double CDFval = CDF(par, j, false);
01431
01432 double CDFlower = CDF(par, j, true);
01433
01434
01435 int qMax = 1;
01436 for (int i = 0; i < nBins; i++) {
01437 if (CDFval > a[i])
01438 qMax = i + 1;
01439 else
01440 break;
01441 }
01442 int qMin = 1;
01443 for (int i = 0; i < nBins; i++) {
01444 if (CDFlower > a[i])
01445 qMin = i + 1;
01446 else
01447 break;
01448 }
01449
01450
01451 if (qMin == qMax) {
01452
01453 if (CDFval < 1)
01454 hist->Fill(CDFval);
01455 else
01456 hist->AddBinContent(qMax);
01457
01458 continue;
01459 }
01460
01461
01462
01463
01464
01465
01466
01467 std::vector<double> prob;
01468
01469 double norm = 1 / double(CDFval - CDFlower);
01470
01471 for (int i = 0; i < (qMax - qMin + 1); i++) {
01472 if (i == 0) {
01473 prob.push_back(norm * (a[qMin] - CDFlower));
01474 continue;
01475 }
01476 if (i == (qMax - qMin)) {
01477 prob.push_back(norm * (CDFval - a[qMax - 1]));
01478 continue;
01479 }
01480
01481 prob.push_back(norm * (a[i] - a[i - 1]));
01482 }
01483
01484 double U = fRandom->Rndm();
01485
01486 for (uint i = 1; i < prob.size(); i++)
01487 prob[i] += prob[i - 1];
01488
01489 for (uint i = 0; i < prob.size(); i++) {
01490
01491 if (U <= prob[i]) {
01492 hist->Fill((a[qMin + i - 1] + a[qMin + i]) / 2.);
01493 break;
01494 }
01495 }
01496 }
01497 }
01498 else {
01499 for (int j = 0; j < n; j++)
01500 hist->Fill(CDF(par, j, false));
01501 }
01502
01503
01504 double chi2 = 0.;
01505 double mk, pk;
01506 double N = double(n);
01507 for (int i = 1; i <= nBins; i++) {
01508 mk = hist->GetBinContent(i);
01509 pk = a[i] - a[i - 1];
01510 chi2 += (mk - N * pk) * (mk - N * pk) / (N * pk);
01511 }
01512
01513 delete hist;
01514
01515 return chi2;
01516 }
01517
01518
01519 double BCModel::GetAvalueFromChi2Johnson(TTree * tree, TH1D * distribution)
01520 {
01521
01522 int nPar = (int)GetNParameters();
01523 std::vector<double> param(nPar);
01524
01525
01526 int nParBranches = -1;
01527 tree->SetBranchAddress("fNParameters", &nParBranches);
01528
01529
01530 tree->GetEntry(0);
01531 if (nParBranches != nPar) {
01532 BCLog::OutError(Form(
01533 "Cannot compute A value: number of parameters in tree (%d)"
01534 "doesn't match number of parameters in model (%d)",
01535 nParBranches, nPar));
01536 return -1.;
01537 }
01538
01539
01540 char * branchName = new char[10 + nPar];
01541
01542
01543 for (int i = 0; i < (int) nPar; i++) {
01544 sprintf(branchName, "fParameter%d", i);
01545 tree->SetBranchAddress(branchName, ¶m[i]);
01546 }
01547
01548
01549 long nEntries = tree->GetEntries();
01550
01551
01552 std::vector<double> randoms(nEntries);
01553 int K = NumberBins();
01554 BCMath::RandomChi2(randoms, K - 1);
01555
01556
01557 int nBigger = 0;
01558 for (int i = 0; i < nEntries; i++) {
01559 tree->GetEntry(i);
01560 double chi2 = GetChi2Johnson(param);
01561 if (distribution != 0)
01562 distribution->Fill(chi2);
01563
01564 if (chi2 > randoms[i])
01565 nBigger++;
01566 }
01567
01568 return nBigger / double(nEntries);
01569 }
01570
01571
01572 double BCModel::GetPvalueFromChi2NDoF(std::vector<double> par, int sigma_index)
01573 {
01574 double ll = LogLikelihood(par);
01575 int n = GetNDataPoints();
01576 int npar = GetNParameters();
01577
01578 double sum_sigma = 0;
01579 for (int i = 0; i < n; i++)
01580 sum_sigma += log(GetDataPoint(i)->GetValue(sigma_index));
01581
01582 double chi2 = -2. * (ll + (double) n / 2. * log(2. * M_PI) + sum_sigma);
01583
01584 fChi2NDoF = chi2 / double(n - npar);
01585 fPValueNDoF = TMath::Prob(chi2, n - npar);
01586
01587 return fPValueNDoF;
01588 }
01589
01590
01591 double BCModel::GetPvalueFromKolmogorov(const std::vector<double>& par,int index)
01592 {
01593 if (flag_discrete) {
01594 BCLog::OutError(Form("BCModel::GetPvalueFromKolmogorov : "
01595 "test defined only for continuous distributions."));
01596 return -1.;
01597 }
01598
01599
01600 std::vector<double> yData = fDataSet->GetDataComponents(index);
01601 TH1D * ECDF = BCMath::ECDF(yData);
01602
01603 int N = GetNDataPoints();
01604
01605
01606 std::set<double> uniqueObservations;
01607 for (int i = 0; i < N; i++)
01608 uniqueObservations.insert(CDF(par, i, false));
01609
01610 int nUnique = uniqueObservations.size();
01611 if (nUnique != ECDF->GetNbinsX() + 1) {
01612 BCLog::OutError(Form("BCModel::GetPvalueFromKolmogorov : "
01613 "Number of unique data points doesn't match (%d vs %d)", nUnique,
01614 ECDF->GetNbinsX() + 1));
01615 return -1.;
01616 }
01617
01618
01619 double distMax = 0.;
01620
01621
01622 double dist = 0.;
01623
01624 std::set<double>::const_iterator iter = uniqueObservations.begin();
01625 for (int iBin = 0; iBin < nUnique; ++iBin) {
01626
01627 dist = TMath::Abs(*iter - ECDF->GetBinContent(iBin + 1));
01628
01629 distMax = TMath::Max(dist, distMax);
01630
01631
01632
01633
01634
01635
01636 ++iter;
01637 }
01638
01639
01640
01641 double z = distMax * sqrt(N);
01642
01643 fPValue = TMath::KolmogorovProb(z);
01644
01645
01646
01647
01648
01649 delete ECDF;
01650
01651 return fPValue;
01652 }
01653
01654
01655 BCH1D * BCModel::CalculatePValue(std::vector<double> par, bool flag_histogram)
01656 {
01657 BCH1D * hist = 0;
01658
01659
01660 BCLog::OutSummary("Do goodness-of-fit-test");
01661
01662
01663 BCGoFTest * goftest = new BCGoFTest("modeltest");
01664
01665
01666 goftest->SetTestModel(this);
01667
01668
01669
01670 if (!goftest->SetTestPoint(par))
01671 return 0;
01672
01673
01674
01675 goftest->MCMCSetFlagFillHistograms(false);
01676
01677
01678 goftest->MCMCSetNChains(fGoFNChains);
01679 goftest->MCMCSetNIterationsMax(fGoFNIterationsMax);
01680 goftest->MCMCSetNIterationsRun(fGoFNIterationsRun);
01681
01682
01683 fPValue = goftest->GetCalculatedPValue(flag_histogram);
01684
01685
01686 if (flag_histogram) {
01687 hist = new BCH1D();
01688 hist->SetHistogram(goftest->GetHistogramLogProb());
01689 }
01690
01691
01692 delete goftest;
01693
01694
01695 return hist;
01696 }
01697
01698
01699 void BCModel::CorrelateDataPointValues(std::vector<double> & )
01700 {
01701
01702 }
01703
01704
01705 double BCModel::HessianMatrixElement(BCParameter * par1, BCParameter * par2, std::vector<double> point)
01706 {
01707
01708 if (point.size() != GetNParameters()) {
01709 BCLog::OutError("BCModel::HessianMatrixElement : Invalid number of entries in the vector.");
01710 return -1;
01711 }
01712
01713
01714 double nsteps = 1e5;
01715 double dx1 = par1->GetRangeWidth() / nsteps;
01716 double dx2 = par2->GetRangeWidth() / nsteps;
01717
01718
01719 std::vector<double> xpp = point;
01720 std::vector<double> xpm = point;
01721 std::vector<double> xmp = point;
01722 std::vector<double> xmm = point;
01723
01724 int idx1 = par1->GetIndex();
01725 int idx2 = par2->GetIndex();
01726
01727 xpp[idx1] += dx1;
01728 xpp[idx2] += dx2;
01729
01730 xpm[idx1] += dx1;
01731 xpm[idx2] -= dx2;
01732
01733 xmp[idx1] -= dx1;
01734 xmp[idx2] += dx2;
01735
01736 xmm[idx1] -= dx1;
01737 xmm[idx2] -= dx2;
01738
01739
01740 double ppp = Likelihood(xpp);
01741 double ppm = Likelihood(xpm);
01742 double pmp = Likelihood(xmp);
01743 double pmm = Likelihood(xmm);
01744
01745
01746 return (ppp + pmm - ppm - pmp) / (4. * dx1 * dx2);
01747 }
01748
01749
01750 void BCModel::FixDataAxis(unsigned int index, bool fixed)
01751 {
01752
01753 if (index > fDataSet->GetDataPoint(0)->GetNValues()) {
01754 BCLog::OutError("BCModel::FixDataAxis : Index out of range.");
01755 return;
01756 }
01757
01758 if (fDataFixedValues.size() == 0)
01759 fDataFixedValues.assign(fDataSet->GetDataPoint(0)->GetNValues(),
01760 false);
01761
01762 fDataFixedValues[index] = fixed;
01763 }
01764
01765
01766 bool BCModel::GetFixedDataAxis(unsigned int index)
01767 {
01768
01769 if (index > fDataSet->GetDataPoint(0)->GetNValues()) {
01770 BCLog::OutError("BCModel::GetFixedDataAxis : Index out of range.");
01771 return false;
01772 }
01773
01774 return fDataFixedValues[index];
01775 }
01776
01777
01778 int BCModel::SetPrior(int index, TF1 * f)
01779 {
01780
01781 if (index < 0 && index >= int(GetNParameters())) {
01782 BCLog::OutError("BCModel::SetPrior : Index out of range.");
01783 return 0;
01784 }
01785
01786 if (fPriorContainer[index])
01787 delete fPriorContainer[index];
01788
01789
01790 fPriorContainer[index] = new TF1(*f);
01791
01792 fPriorContainerConstant[index] = false;
01793
01794 RecalculatePriorConstant();
01795
01796
01797 ResetResults();
01798
01799
01800 return 1;
01801 }
01802
01803
01804 int BCModel::SetPrior(const char * name, TF1 * f)
01805 {
01806
01807 int index = -1;
01808 for (unsigned int i = 0; i < GetNParameters(); i++)
01809 if (name == GetParameter(i)->GetName())
01810 index = i;
01811
01812
01813 return SetPrior(index, f);
01814 }
01815
01816
01817 int BCModel::SetPriorDelta(int index, double value)
01818 {
01819
01820 SetParameterRange(index, value, value);
01821
01822
01823 return SetPriorConstant(index);
01824 }
01825
01826
01827 int BCModel::SetPriorDelta(const char* name, double value)
01828 {
01829
01830 int index = -1;
01831 for (unsigned int i = 0; i < GetNParameters(); i++)
01832 if (name == GetParameter(i)->GetName())
01833 index = i;
01834
01835
01836 return SetPriorDelta(index, value);
01837 }
01838
01839
01840 int BCModel::SetPriorGauss(int index, double mean, double sigma)
01841 {
01842
01843 if (index < 0 || index >= int(GetNParameters())) {
01844 BCLog::OutError("BCModel::SetPriorGauss : Index out of range.");
01845 return 0;
01846 }
01847
01848
01849 TF1 * f = new TF1(Form("prior_%s", GetParameter(index)->GetName().c_str()),
01850 "1./sqrt(2.*TMath::Pi())/[1] * exp(- (x-[0])*(x-[0])/2./[1]/[1])",
01851 GetParameter(index)->GetLowerLimit(),
01852 GetParameter(index)->GetUpperLimit());
01853 f->SetParameter(0, mean);
01854 f->SetParameter(1, sigma);
01855
01856
01857 return SetPrior(index, f);
01858 }
01859
01860
01861 int BCModel::SetPriorGauss(const char* name, double mean, double sigma)
01862 {
01863
01864 int index = -1;
01865 for (unsigned int i = 0; i < GetNParameters(); i++)
01866 if (name == GetParameter(i)->GetName())
01867 index = i;
01868
01869
01870 return SetPriorGauss(index, mean, sigma);
01871 }
01872
01873
01874 int BCModel::SetPriorGauss(int index, double mean, double sigmadown, double sigmaup)
01875 {
01876
01877 if (index < 0 || index >= int(GetNParameters())) {
01878 BCLog::OutError("BCModel::SetPriorGauss : Index out of range.");
01879 return 0;
01880 }
01881
01882
01883 TF1 * f = new TF1(Form("prior_%s", GetParameter(index)->GetName().c_str()),
01884 BCMath::SplitGaussian,
01885 GetParameter(index)->GetLowerLimit(),
01886 GetParameter(index)->GetUpperLimit(),
01887 3);
01888 f->SetParameter(0, mean);
01889 f->SetParameter(1, sigmadown);
01890 f->SetParameter(2, sigmaup);
01891
01892
01893 return SetPrior(index, f);
01894
01895 return 0;
01896 }
01897
01898
01899 int BCModel::SetPriorGauss(const char * name, double mean, double sigmadown, double sigmaup)
01900 {
01901
01902 int index = -1;
01903 for (unsigned int i = 0; i < GetNParameters(); i++)
01904 if (name == GetParameter(i)->GetName())
01905 index = i;
01906
01907
01908 return SetPriorGauss(index, mean, sigmadown, sigmaup);
01909 }
01910
01911
01912 int BCModel::SetPrior(int index, TH1 * h, bool interpolate)
01913 {
01914
01915 if (index < 0 && index >= int(GetNParameters())) {
01916 BCLog::OutError("BCModel::SetPrior : Index out of range.");
01917 return 0;
01918 }
01919
01920
01921 if(h) {
01922
01923
01924 if (h->GetDimension() != 1) {
01925 BCLog::OutError(Form("BCModel::SetPrior : Histogram given for parameter %d is not 1D.",index));
01926 return 0;
01927 }
01928
01929
01930 h->Scale(1./h->Integral("width"));
01931
01932 if(fPriorContainer[index])
01933 delete fPriorContainer[index];
01934
01935
01936 fPriorContainer[index] = new TH1(*h);
01937
01938 if (interpolate)
01939 fPriorContainerInterpolate[index] = true;
01940
01941 fPriorContainerConstant[index] = false;
01942 }
01943
01944 RecalculatePriorConstant();
01945
01946
01947 ResetResults();
01948
01949
01950 return 1;
01951 }
01952
01953
01954 int BCModel::SetPrior(const char * name, TH1 * h, bool interpolate)
01955 {
01956
01957 int index = -1;
01958 for (unsigned int i = 0; i < GetNParameters(); i++)
01959 if (name == GetParameter(i)->GetName())
01960 index = i;
01961
01962
01963 return SetPrior(index, h, interpolate);
01964 }
01965
01966
01967 int BCModel::SetPriorConstant(int index)
01968 {
01969
01970 if (index < 0 && index >= int(GetNParameters())) {
01971 BCLog::OutError("BCModel::SetPriorConstant : Index out of range.");
01972 return 0;
01973 }
01974
01975 if(fPriorContainer[index]) {
01976 delete fPriorContainer[index];
01977 fPriorContainer[index] = 0;
01978 }
01979
01980
01981 fPriorContainerConstant[index] = true;
01982
01983 RecalculatePriorConstant();
01984
01985
01986 ResetResults();
01987
01988
01989 return 1;
01990 }
01991
01992
01993 int BCModel::SetPriorConstant(const char* name)
01994 {
01995
01996 int index = -1;
01997 for (unsigned int i = 0; i < GetNParameters(); i++) {
01998 if (name == GetParameter(i)->GetName()) {
01999 index = i;
02000 break;
02001 }
02002 }
02003
02004 if (index == -1) {
02005 BCLog::OutError(Form(
02006 "BCModel::SetPriorConstant : parameter '%s' doesn't exist.", name));
02007 return 0;
02008 }
02009
02010 return SetPriorConstant(index);
02011 }
02012
02013
02014 int BCModel::SetPriorConstantAll()
02015 {
02016
02017 int nPar = GetNParameters();
02018
02019 if (nPar == 0)
02020 BCLog::OutWarning("BCModel::SetPriorConstantAll : No parameters defined.");
02021
02022
02023 for (int i = 0; i < nPar; ++i) {
02024 if (fPriorContainer[i]) {
02025 delete fPriorContainer[i];
02026 fPriorContainer[i]=0;
02027 }
02028 fPriorContainerConstant[i] = true;
02029 }
02030
02031 RecalculatePriorConstant();
02032
02033
02034 ResetResults();
02035
02036
02037 return 1;
02038 }
02039
02040
02041 void BCModel::RecalculatePriorConstant()
02042 {
02043 fPriorConstantValue = 0.;
02044
02045
02046 int npar = GetNParameters();
02047
02048 int nconstant = 0;
02049
02050 for (int i=0; i<npar; ++i)
02051 if (fPriorContainerConstant[i]) {
02052
02053 if (GetParameter(i)->GetRangeWidth() > 0)
02054 fPriorConstantValue -= log(GetParameter(i)->GetRangeWidth());
02055
02056 ++nconstant;
02057 }
02058
02059 if (nconstant == npar)
02060 fPriorConstantAll = true;
02061 else
02062 fPriorConstantAll = false;
02063 }
02064
02065
02066 int BCModel::SetParameterRange(int index, double parmin, double parmax)
02067 {
02068
02069 if (index < 0 || index >= int(GetNParameters())) {
02070 BCLog::OutError("BCModel::SetParameterRange : Index out of range.");
02071 return 0;
02072 }
02073
02074
02075 GetParameter(index)->SetLowerLimit(parmin);
02076 GetParameter(index)->SetUpperLimit(parmax);
02077 fMCMCBoundaryMin[index] = parmin;
02078 fMCMCBoundaryMax[index] = parmax;
02079
02080
02081 ResetResults();
02082
02083
02084 return 1;
02085 }
02086
02087
02088 int BCModel::ResetResults()
02089 {
02090 BCIntegrate::IntegrateResetResults();
02091
02092 BCEngineMCMC::MCMCResetResults();
02093
02094
02095 return 1;
02096 }
02097
02098
02099 void BCModel::PrintSummary()
02100 {
02101 int nparameters = GetNParameters();
02102
02103
02104 std::cout << std::endl
02105 << " ---------------------------------" << std::endl
02106 << " Model : " << fName.data() << std::endl
02107 << " ---------------------------------" << std::endl
02108 << " Index : " << fIndex << std::endl
02109 << " Number of parameters : " << nparameters << std::endl << std::endl
02110 << " - Parameters : " << std::endl << std::endl;
02111
02112
02113 for (int i = 0; i < nparameters; i++)
02114 fParameterSet->at(i)->PrintSummary();
02115
02116
02117 if (GetBestFitParameters().size() > 0) {
02118 std::cout << std::endl << " - Best fit parameters :" << std::endl << std::endl;
02119
02120 for (int i = 0; i < nparameters; i++) {
02121 std::cout << " " << fParameterSet->at(i)->GetName().data()
02122 << " = " << GetBestFitParameter(i) << " (overall)" << std::endl;
02123 if ((int) fBestFitParametersMarginalized.size() == nparameters)
02124 std::cout << " "
02125 << fParameterSet->at(i)->GetName().data() << " = "
02126 << GetBestFitParameterMarginalized(i)
02127 << " (marginalized)" << std::endl;
02128 }
02129 }
02130
02131 std::cout << std::endl;
02132
02133
02134 if (fPValue >= 0) {
02135 double likelihood = Likelihood(GetBestFitParameters());
02136 std::cout << " - Model testing:" << std::endl << std::endl
02137 << " p(data|lambda*) = " << likelihood << std::endl
02138 << " p-value = " << fPValue << std::endl << std::endl;
02139 }
02140
02141
02142 if (fNormalization > 0)
02143 std::cout << " Normalization : " << fNormalization << std::endl;
02144 }
02145
02146
02147 void BCModel::PrintResults(const char * file)
02148 {
02149
02150
02151
02152 ofstream ofi(file);
02153
02154
02155 if (!ofi.is_open()) {
02156 std::cerr << "Couldn't open file " << file << std::endl;
02157 return;
02158 }
02159
02160
02161 int npar = MCMCGetNParameters();
02162 int nchains = MCMCGetNChains();
02163
02164
02165 bool flag_conv = MCMCGetNIterationsConvergenceGlobal() > 0;
02166
02167 ofi << std::endl
02168 << " -----------------------------------------------------" << std::endl
02169 << " Summary" << std::endl
02170 << " -----------------------------------------------------" << std::endl
02171 << std::endl;
02172
02173 ofi << " Model summary" << std::endl << " =============" << std::endl
02174 << " Model: " << fName.data() << std::endl
02175 << " Number of parameters: " << npar << std::endl
02176 << " List of Parameters and ranges:" << std::endl;
02177 for (int i = 0; i < npar; ++i)
02178 ofi << " (" << i << ") Parameter \""
02179 << fParameterSet->at(i)->GetName().data() << "\"" << ": "
02180 << "(" << fParameterSet->at(i)->GetLowerLimit() << ", "
02181 << fParameterSet->at(i)->GetUpperLimit() << ")" << std::endl;
02182 ofi << std::endl;
02183
02184 ofi << " Results of the optimization" << std::endl
02185 << " ===========================" << std::endl
02186 << " Optimization algorithm used: "
02187 << DumpUsedOptimizationMethod()<< std::endl;
02188
02189 if (int(fBestFitParameters.size()) > 0) {
02190 ofi << " List of parameters and global mode:" << std::endl;
02191 for (int i = 0; i < npar; ++i) {
02192 ofi << " (" << i << ") Parameter \""
02193 << fParameterSet->at(i)->GetName().data() << "\": "
02194 << fBestFitParameters[i];
02195 if (int(fBestFitParameterErrors.size()) == npar)
02196 if(fBestFitParameterErrors[i]>=0.)
02197 ofi << " +- " << fBestFitParameterErrors[i];
02198 ofi << std::endl;
02199 }
02200 ofi << std::endl;
02201 }
02202 else {
02203 ofi << " No best fit information available." << std::endl;
02204 ofi << std::endl;
02205 }
02206
02207 if (fPValue >= 0.) {
02208 ofi << " Results of the model test" << std::endl
02209 << " =========================" << std::endl
02210 << " p-value at global mode: " << fPValue << std::endl << std::endl;
02211 }
02212
02213 if (fNormalization >= 0.) {
02214 ofi << " Results of the normalization" << std::endl
02215 << " ============================" << std::endl
02216 << " Integration method used:"
02217 << DumpIntegrationMethod() << std::endl;
02218 ofi << " Normalization factor: " << fNormalization << std::endl << std::endl;
02219 }
02220
02221
02222 if (!flag_conv && fMCMCFlagRun)
02223 ofi << " WARNING: the Markov Chain did not converge!" << std::endl
02224 << " Be cautious using the following results!" << std::endl
02225 << std::endl;
02226
02227
02228 if (fMCMCFlagRun && fMCMCFlagFillHistograms) {
02229 ofi << " Results of the marginalization" << std::endl
02230 << " ==============================" << std::endl
02231 << " List of parameters and properties of the marginalized"
02232 << std::endl << " distributions:" << std::endl;
02233 for (int i = 0; i < npar; ++i) {
02234 if (!fMCMCFlagsFillHistograms[i])
02235 continue;
02236
02237 BCH1D * bch1d = GetMarginalized(fParameterSet->at(i));
02238
02239 ofi << " (" << i << ") Parameter \""
02240 << fParameterSet->at(i)->GetName().data() << "\":" << std::endl
02241
02242 << " Mean +- sqrt(V): " << std::setprecision(4)
02243 << bch1d->GetMean() << " +- " << std::setprecision(4)
02244 << bch1d->GetRMS() << std::endl
02245
02246 << " Median +- central 68% interval: "
02247 << std::setprecision(4) << bch1d->GetMedian() << " + "
02248 << std::setprecision(4) << bch1d->GetQuantile(0.84) - bch1d->GetMedian()
02249 << " - " << std::setprecision(4)
02250 << bch1d->GetMedian() - bch1d->GetQuantile(0.16) << std::endl
02251
02252 << " (Marginalized) mode: " << bch1d->GetMode() << std::endl;
02253
02254 ofi << " 5% quantile: " << std::setprecision(4)
02255 << bch1d->GetQuantile(0.05) << std::endl
02256 << " 10% quantile: " << std::setprecision(4)
02257 << bch1d->GetQuantile(0.10) << std::endl
02258 << " 16% quantile: " << std::setprecision(4)
02259 << bch1d->GetQuantile(0.16) << std::endl
02260 << " 84% quantile: " << std::setprecision(4)
02261 << bch1d->GetQuantile(0.85) << std::endl
02262 << " 90% quantile: " << std::setprecision(4)
02263 << bch1d->GetQuantile(0.90) << std::endl
02264 << " 95% quantile: " << std::setprecision(4)
02265 << bch1d->GetQuantile(0.95) << std::endl;
02266
02267 std::vector<double> v;
02268 v = bch1d->GetSmallestIntervals(0.68);
02269 int ninter = int(v.size());
02270 ofi << " Smallest interval(s) containing 68% and local modes:"
02271 << std::endl;
02272 for (int j = 0; j < ninter; j += 5)
02273 ofi << " (" << v[j] << ", " << v[j + 1]
02274 << ") (local mode at " << v[j + 3] << " with rel. height "
02275 << v[j + 2] << "; rel. area " << v[j + 4] << ")"
02276 << std::endl;
02277 ofi << std::endl;
02278 }
02279 }
02280 if (fMCMCFlagRun) {
02281 ofi << " Status of the MCMC" << std::endl << " ==================" << std::endl
02282 << " Convergence reached: " << (flag_conv ? "yes" : "no")
02283 << std::endl;
02284
02285 if (flag_conv)
02286 ofi << " Number of iterations until convergence: "
02287 << MCMCGetNIterationsConvergenceGlobal() << std::endl;
02288 else
02289 ofi << " WARNING: the Markov Chain did not converge! Be\n"
02290 << " cautious using the following results!" << std::endl
02291 << std::endl;
02292 ofi << " Number of chains: " << MCMCGetNChains() << std::endl
02293 << " Number of iterations per chain: " << MCMCGetNIterationsRun() << std::endl
02294 << " Average efficiencies:" << std::endl;
02295
02296 std::vector<double> efficiencies;
02297 efficiencies.assign(npar, 0.);
02298
02299 for (int ipar = 0; ipar < npar; ++ipar)
02300 for (int ichain = 0; ichain < nchains; ++ichain) {
02301 int index = ichain * npar + ipar;
02302 efficiencies[ipar] +=
02303 double(MCMCGetNTrialsTrue().at(index)) / double(MCMCGetNTrialsTrue().at(index)
02304 + MCMCGetNTrialsFalse().at(index)) / double(nchains) * 100.;
02305 }
02306
02307 for (int ipar = 0; ipar < npar; ++ipar)
02308 ofi << " (" << ipar << ") Parameter \""
02309 << fParameterSet->at(ipar)->GetName().data() << "\": "
02310 << efficiencies.at(ipar) << "%" << std::endl;
02311 ofi << std::endl;
02312 }
02313
02314 ofi << " -----------------------------------------------------" << std::endl
02315 << " Notation:" << std::endl
02316 << " Mean : mean value of the marg. pdf" << std::endl
02317 << " Median : maximum of the marg. pdf" << std::endl
02318 << " Marg. mode : most probable value of the marg. pdf" << std::endl
02319 << " V : Variance of the marg. pdf" << std::endl
02320 << " Quantiles : most commonly used quantiles" <<std::endl
02321 << " -----------------------------------------------------" << std::endl
02322 << std::endl;
02323
02324
02325
02326 }
02327
02328
02329 void BCModel::PrintShortFitSummary(int chi2flag)
02330 {
02331 BCLog::OutSummary("---------------------------------------------------");
02332 BCLog::OutSummary(Form("Fit summary for model \'%s\':", GetName().data()));
02333 BCLog::OutSummary(Form(" Number of parameters: Npar = %i", GetNParameters()));
02334 BCLog::OutSummary(Form(" Number of data points: Ndata = %i", GetNDataPoints()));
02335 BCLog::OutSummary(" Number of degrees of freedom:");
02336 BCLog::OutSummary(Form(" NDoF = Ndata - Npar = %i", GetNDataPoints() - GetNParameters()));
02337
02338 BCLog::OutSummary(" Best fit parameters (global):");
02339 for (unsigned int i = 0; i < GetNParameters(); ++i)
02340 BCLog::OutSummary(Form(" %s = %.3g", GetParameter(i)->GetName().data(), GetBestFitParameter(i)));
02341
02342 BCLog::OutSummary(" Goodness-of-fit test:");
02343 BCLog::OutSummary(Form(" p-value = %.3g", GetPValue()));
02344 if (chi2flag) {
02345 BCLog::OutSummary(Form(" p-value corrected for NDoF = %.3g", GetPValueNDoF()));
02346 BCLog::OutSummary(Form(" chi2 / NDoF = %.3g", GetChi2NDoF()));
02347 }
02348 BCLog::OutSummary("---------------------------------------------------");
02349 }
02350
02351
02352 void BCModel::PrintHessianMatrix(std::vector<double> parameters)
02353 {
02354
02355 if (parameters.size() != GetNParameters()) {
02356 BCLog::OutError("BCModel::PrintHessianMatrix : Invalid number of entries in the vector");
02357 return;
02358 }
02359
02360
02361 std::cout << std::endl << " Hessian matrix elements: " << std::endl << " Point: ";
02362
02363 for (int i = 0; i < int(parameters.size()); i++)
02364 std::cout << parameters.at(i) << " ";
02365 std::cout << std::endl;
02366
02367
02368 for (unsigned int i = 0; i < GetNParameters(); i++)
02369 for (unsigned int j = 0; j < i; j++) {
02370
02371 double hessianmatrixelement = HessianMatrixElement(
02372 fParameterSet->at(i), fParameterSet->at(j), parameters);
02373
02374
02375 std::cout << " " << i << " " << j << " : " << hessianmatrixelement << std::endl;
02376 }
02377 }
02378
02379
02380 BCDataPoint * BCModel::VectorToDataPoint(const std::vector<double> &data)
02381 {
02382 int sizeofvector = int(data.size());
02383 BCDataPoint * datapoint = new BCDataPoint(sizeofvector);
02384 datapoint->SetValues(data);
02385 return datapoint;
02386 }
02387
02388
02389 int BCModel::CompareStrings(const char * string1, const char * string2)
02390 {
02391 int flag_same = 0;
02392
02393 if (strlen(string1) != strlen(string2))
02394 return -1;
02395
02396 for (int i = 0; i < int(strlen(string1)); i++)
02397 if (string1[i] != string2[i])
02398 flag_same = -1;
02399
02400 return flag_same;
02401 }
02402
02403 void BCModel::StoreMode()
02404 {
02405
02406 double probmax = fMCMCprobMax.at(0);
02407 int probmaxindex = 0;
02408
02409
02410 for (int i = 1; i < fMCMCNChains; ++i) {
02411 if (fMCMCprobMax.at(i) > probmax) {
02412 probmax = fMCMCprobMax.at(i);
02413 probmaxindex = i;
02414 }
02415 }
02416
02417
02418 if (fBestFitParameters.empty() || probmax > LogEval(fBestFitParameters)) {
02419 fBestFitParameters.assign(fMCMCNParameters, 0.0);
02420 for (int i = 0; i < fMCMCNParameters; ++i)
02421 fBestFitParameters[i] = fMCMCxMax[probmaxindex * fMCMCNParameters + i];
02422
02423 SetOptimizationMethodMode(BCIntegrate::kOptMetropolis);
02424 }
02425 }
02426
02427
02428