00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064 #include "RooDataHist.h"
00065 #include "RooDataSet.h"
00066 #include "RooGlobalFunc.h"
00067 #include "RooNLLVar.h"
00068 #include "RooRealVar.h"
00069 #include "RooAbsData.h"
00070 #include "RooWorkspace.h"
00071
00072 #include "TH1.h"
00073
00074 #include "RooStats/HybridCalculatorOriginal.h"
00075
00076 ClassImp(RooStats::HybridCalculatorOriginal)
00077
00078 using namespace RooStats;
00079
00080
00081
00082 HybridCalculatorOriginal::HybridCalculatorOriginal(const char *name) :
00083 TNamed(name,name),
00084 fSbModel(0),
00085 fBModel(0),
00086 fObservables(0),
00087 fNuisanceParameters(0),
00088 fPriorPdf(0),
00089 fData(0),
00090 fGenerateBinned(false),
00091 fUsePriorPdf(false), fTmpDoExtended(true)
00092 {
00093
00094
00095 SetTestStatistic(1);
00096 SetNumberOfToys(1000);
00097 }
00098
00099
00100
00101 HybridCalculatorOriginal::HybridCalculatorOriginal( RooAbsPdf& sbModel,
00102 RooAbsPdf& bModel,
00103 RooArgList& observables,
00104 const RooArgSet* nuisance_parameters,
00105 RooAbsPdf* priorPdf ,
00106 bool GenerateBinned,
00107 int testStatistics,
00108 int numToys) :
00109 fSbModel(&sbModel),
00110 fBModel(&bModel),
00111 fNuisanceParameters(nuisance_parameters),
00112 fPriorPdf(priorPdf),
00113 fData(0),
00114 fGenerateBinned(GenerateBinned),
00115 fUsePriorPdf(false),
00116 fTmpDoExtended(true)
00117 {
00118
00119
00120
00121
00122
00123
00124 fObservables = new RooArgList(observables);
00125
00126
00127
00128
00129
00130
00131 SetTestStatistic(testStatistics);
00132 SetNumberOfToys(numToys);
00133
00134 if (priorPdf) UseNuisance(true);
00135
00136
00137
00138 }
00139
00140
00141 HybridCalculatorOriginal::HybridCalculatorOriginal( RooAbsData & data,
00142 RooAbsPdf& sbModel,
00143 RooAbsPdf& bModel,
00144 const RooArgSet* nuisance_parameters,
00145 RooAbsPdf* priorPdf,
00146 bool GenerateBinned,
00147 int testStatistics,
00148 int numToys) :
00149 fSbModel(&sbModel),
00150 fBModel(&bModel),
00151 fObservables(0),
00152 fNuisanceParameters(nuisance_parameters),
00153 fPriorPdf(priorPdf),
00154 fData(&data),
00155 fGenerateBinned(GenerateBinned),
00156 fUsePriorPdf(false),
00157 fTmpDoExtended(true)
00158 {
00159
00160
00161
00162
00163
00164
00165 SetTestStatistic(testStatistics);
00166 SetNumberOfToys(numToys);
00167
00168 if (priorPdf) UseNuisance(true);
00169 }
00170
00171
00172
00173 HybridCalculatorOriginal::HybridCalculatorOriginal( RooAbsData& data,
00174 const ModelConfig& sbModel,
00175 const ModelConfig& bModel,
00176 bool GenerateBinned,
00177 int testStatistics,
00178 int numToys) :
00179 fSbModel(sbModel.GetPdf()),
00180 fBModel(bModel.GetPdf()),
00181 fObservables(0),
00182 fNuisanceParameters((sbModel.GetNuisanceParameters()) ? sbModel.GetNuisanceParameters() : bModel.GetNuisanceParameters()),
00183 fPriorPdf((sbModel.GetPriorPdf()) ? sbModel.GetPriorPdf() : bModel.GetPriorPdf()),
00184 fData(&data),
00185 fGenerateBinned(GenerateBinned),
00186 fUsePriorPdf(false),
00187 fTmpDoExtended(true)
00188 {
00189
00190
00191
00192
00193
00194 if (fPriorPdf) UseNuisance(true);
00195
00196 SetTestStatistic(testStatistics);
00197 SetNumberOfToys(numToys);
00198 }
00199
00200
00201
00202 HybridCalculatorOriginal::~HybridCalculatorOriginal()
00203 {
00204
00205 if (fObservables) delete fObservables;
00206 }
00207
00208
00209
00210 void HybridCalculatorOriginal::SetNullModel(const ModelConfig& model)
00211 {
00212
00213 fBModel = model.GetPdf();
00214
00215 if (!fPriorPdf) fPriorPdf = model.GetPriorPdf();
00216 if (!fNuisanceParameters) fNuisanceParameters = model.GetNuisanceParameters();
00217 }
00218
00219 void HybridCalculatorOriginal::SetAlternateModel(const ModelConfig& model)
00220 {
00221
00222 fSbModel = model.GetPdf();
00223 fPriorPdf = model.GetPriorPdf();
00224 fNuisanceParameters = model.GetNuisanceParameters();
00225 }
00226
00227 void HybridCalculatorOriginal::SetTestStatistic(int index)
00228 {
00229
00230
00231
00232
00233
00234 fTestStatisticsIdx = index;
00235 }
00236
00237
00238
00239 HybridResult* HybridCalculatorOriginal::Calculate(TH1& data, unsigned int nToys, bool usePriors) const
00240 {
00241
00242
00243
00244 TString dataHistName = GetName(); dataHistName += "_roodatahist";
00245 RooDataHist dataHist(dataHistName,"Data distribution as RooDataHist converted from TH1",*fObservables,&data);
00246
00247 HybridResult* result = Calculate(dataHist,nToys,usePriors);
00248
00249 return result;
00250 }
00251
00252
00253
00254 HybridResult* HybridCalculatorOriginal::Calculate(RooAbsData& data, unsigned int nToys, bool usePriors) const
00255 {
00256
00257
00258 double testStatData = 0;
00259 if ( fTestStatisticsIdx==2 ) {
00260
00261 double nEvents = data.sumEntries();
00262 testStatData = nEvents;
00263 } else if ( fTestStatisticsIdx==3 ) {
00264
00265 if ( fTmpDoExtended ) {
00266 RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,data,RooFit::CloneData(false),RooFit::Extended());
00267 fSbModel->fitTo(data,RooFit::Hesse(false),RooFit::Strategy(0),RooFit::Extended());
00268 double sb_nll_val = sb_nll.getVal();
00269 RooNLLVar b_nll("b_nll","b_nll",*fBModel,data,RooFit::CloneData(false),RooFit::Extended());
00270 fBModel->fitTo(data,RooFit::Hesse(false),RooFit::Strategy(0),RooFit::Extended());
00271 double b_nll_val = b_nll.getVal();
00272 double m2lnQ = 2*(sb_nll_val-b_nll_val);
00273 testStatData = m2lnQ;
00274 } else {
00275 RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,data,RooFit::CloneData(false));
00276 fSbModel->fitTo(data,RooFit::Hesse(false),RooFit::Strategy(0));
00277 double sb_nll_val = sb_nll.getVal();
00278 RooNLLVar b_nll("b_nll","b_nll",*fBModel,data,RooFit::CloneData(false));
00279 fBModel->fitTo(data,RooFit::Hesse(false),RooFit::Strategy(0));
00280 double b_nll_val = b_nll.getVal();
00281 double m2lnQ = 2*(sb_nll_val-b_nll_val);
00282 testStatData = m2lnQ;
00283 }
00284 } else if ( fTestStatisticsIdx==1 ) {
00285
00286 if ( fTmpDoExtended ) {
00287 RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,data,RooFit::Extended());
00288 RooNLLVar b_nll("b_nll","b_nll",*fBModel,data,RooFit::Extended());
00289 double m2lnQ = 2*(sb_nll.getVal()-b_nll.getVal());
00290 testStatData = m2lnQ;
00291 } else {
00292 RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,data);
00293 RooNLLVar b_nll("b_nll","b_nll",*fBModel,data);
00294 double m2lnQ = 2*(sb_nll.getVal()-b_nll.getVal());
00295 testStatData = m2lnQ;
00296 }
00297 }
00298
00299 std::cout << "Test statistics has been evaluated for data\n";
00300
00301 HybridResult* result = Calculate(nToys,usePriors);
00302 result->SetDataTestStatistics(testStatData);
00303
00304 return result;
00305 }
00306
00307
00308
00309 HybridResult* HybridCalculatorOriginal::Calculate(unsigned int nToys, bool usePriors) const
00310 {
00311 std::vector<double> bVals;
00312 bVals.reserve(nToys);
00313
00314 std::vector<double> sbVals;
00315 sbVals.reserve(nToys);
00316
00317 RunToys(bVals,sbVals,nToys,usePriors);
00318
00319 HybridResult* result;
00320
00321 TString name = "HybridResult_" + TString(GetName() );
00322
00323 if ( fTestStatisticsIdx==2 )
00324 result = new HybridResult(name,sbVals,bVals,false);
00325 else
00326 result = new HybridResult(name,sbVals,bVals);
00327
00328 return result;
00329 }
00330
00331
00332
00333 void HybridCalculatorOriginal::RunToys(std::vector<double>& bVals, std::vector<double>& sbVals, unsigned int nToys, bool usePriors) const
00334 {
00335
00336 std::cout << "HybridCalculatorOriginal: run " << nToys << " toy-MC experiments\n";
00337 std::cout << "with test statistics index: " << fTestStatisticsIdx << "\n";
00338 if (usePriors) std::cout << "marginalize nuisance parameters \n";
00339
00340 assert(nToys > 0);
00341 assert(fBModel);
00342 assert(fSbModel);
00343 if (usePriors) {
00344 assert(fPriorPdf);
00345 assert(fNuisanceParameters);
00346 }
00347
00348 std::vector<double> parameterValues;
00349
00350 int nParameters = (fNuisanceParameters) ? fNuisanceParameters->getSize() : 0;
00351 RooArgList parametersList("parametersList");
00352 if (usePriors && nParameters>0) {
00353 parametersList.add(*fNuisanceParameters);
00354 parameterValues.resize(nParameters);
00355 for (int iParameter=0; iParameter<nParameters; iParameter++) {
00356 RooRealVar* oneParam = (RooRealVar*) parametersList.at(iParameter);
00357 parameterValues[iParameter] = oneParam->getVal();
00358 }
00359 }
00360
00361
00362
00363 RooArgSet originalSbParams;
00364 RooArgSet originalBParams;
00365 if (fTestStatisticsIdx == 3) {
00366 RooArgSet * sbparams = fSbModel->getParameters(*fObservables);
00367 RooArgSet * bparams = fBModel->getParameters(*fObservables);
00368 if (sbparams) originalSbParams.addClone(*sbparams);
00369 if (bparams) originalBParams.addClone(*bparams);
00370 delete sbparams;
00371 delete bparams;
00372
00373
00374 }
00375
00376
00377 for (unsigned int iToy=0; iToy<nToys; iToy++) {
00378
00379
00380
00381 if ( iToy%500==0 ) {
00382 std::cout << "....... toy number " << iToy << " / " << nToys << std::endl;
00383 }
00384
00385
00386 if (usePriors && nParameters>0) {
00387
00388 RooDataSet* tmpValues = (RooDataSet*) fPriorPdf->generate(*fNuisanceParameters,1);
00389 for (int iParameter=0; iParameter<nParameters; iParameter++) {
00390 RooRealVar* oneParam = (RooRealVar*) parametersList.at(iParameter);
00391 oneParam->setVal(tmpValues->get()->getRealValue(oneParam->GetName()));
00392 }
00393 delete tmpValues;
00394 }
00395
00396
00397
00398 RooAbsData* bData;
00399 if (fGenerateBinned)
00400 bData = static_cast<RooAbsData*> (fBModel->generateBinned(*fObservables,RooFit::Extended()));
00401 else {
00402 if ( fTmpDoExtended ) bData = static_cast<RooAbsData*> (fBModel->generate(*fObservables,RooFit::Extended()));
00403 else bData = static_cast<RooAbsData*> (fBModel->generate(*fObservables,1));
00404 }
00405
00406
00407 bool bIsEmpty = false;
00408 if (bData==NULL) {
00409 bIsEmpty = true;
00410
00411 RooDataSet* bDataDummy=new RooDataSet("bDataDummy","empty dataset",*fObservables);
00412 bData = static_cast<RooAbsData*>(new RooDataHist ("bDataEmpty","",*fObservables,*bDataDummy));
00413 delete bDataDummy;
00414 }
00415
00416
00417 RooAbsData* sbData;
00418 if (fGenerateBinned)
00419 sbData = static_cast<RooAbsData*> (fSbModel->generateBinned(*fObservables,RooFit::Extended()));
00420 else {
00421 if ( fTmpDoExtended ) sbData = static_cast<RooAbsData*> (fSbModel->generate(*fObservables,RooFit::Extended()));
00422 else sbData = static_cast<RooAbsData*> (fSbModel->generate(*fObservables,1));
00423 }
00424
00425
00426 bool sbIsEmpty = false;
00427 if (sbData==NULL) {
00428 sbIsEmpty = true;
00429
00430 RooDataSet* sbDataDummy=new RooDataSet("sbDataDummy","empty dataset",*fObservables);
00431 sbData = static_cast<RooAbsData*>(new RooDataHist ("sbDataEmpty","",*fObservables,*sbDataDummy));
00432 delete sbDataDummy;
00433 }
00434
00435
00436 if (usePriors && nParameters>0) {
00437 for (int iParameter=0; iParameter<nParameters; iParameter++) {
00438 RooRealVar* oneParam = (RooRealVar*) parametersList.at(iParameter);
00439 oneParam->setVal(parameterValues[iParameter]);
00440 }
00441 }
00442
00443
00444 for (int hypoTested=0; hypoTested<=1; hypoTested++) {
00445 RooAbsData* dataToTest = sbData;
00446 bool dataIsEmpty = sbIsEmpty;
00447 if ( hypoTested==1 ) { dataToTest = bData; dataIsEmpty = bIsEmpty; }
00448
00449 if ( fTestStatisticsIdx==2 ) {
00450 double nEvents = 0;
00451 if ( !dataIsEmpty ) nEvents = dataToTest->numEntries();
00452 if ( hypoTested==0 ) sbVals.push_back(nEvents);
00453 else bVals.push_back(nEvents);
00454 } else if ( fTestStatisticsIdx==3 ) {
00455 if ( fTmpDoExtended ) {
00456 RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,*dataToTest,RooFit::CloneData(false),RooFit::Extended());
00457 fSbModel->fitTo(*dataToTest,RooFit::PrintLevel(-1), RooFit::Hesse(false),RooFit::Strategy(0),RooFit::Extended());
00458 double sb_nll_val = sb_nll.getVal();
00459 RooNLLVar b_nll("b_nll","b_nll",*fBModel,*dataToTest,RooFit::CloneData(false),RooFit::Extended());
00460 fBModel->fitTo(*dataToTest,RooFit::PrintLevel(-1),RooFit::Hesse(false),RooFit::Strategy(0),RooFit::Extended());
00461 double b_nll_val = b_nll.getVal();
00462 double m2lnQ = -2*(b_nll_val-sb_nll_val);
00463 if ( hypoTested==0 ) sbVals.push_back(m2lnQ);
00464 else bVals.push_back(m2lnQ);
00465 } else {
00466 RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,*dataToTest,RooFit::CloneData(false));
00467 fSbModel->fitTo(*dataToTest,RooFit::PrintLevel(-1), RooFit::Hesse(false),RooFit::Strategy(0));
00468 double sb_nll_val = sb_nll.getVal();
00469 RooNLLVar b_nll("b_nll","b_nll",*fBModel,*dataToTest,RooFit::CloneData(false));
00470 fBModel->fitTo(*dataToTest,RooFit::PrintLevel(-1), RooFit::Hesse(false),RooFit::Strategy(0));
00471 double b_nll_val = b_nll.getVal();
00472 double m2lnQ = -2*(b_nll_val-sb_nll_val);
00473 if ( hypoTested==0 ) sbVals.push_back(m2lnQ);
00474 else bVals.push_back(m2lnQ);
00475 }
00476 } else if ( fTestStatisticsIdx==1 ) {
00477 if ( fTmpDoExtended ) {
00478 RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,*dataToTest,RooFit::CloneData(false),RooFit::Extended());
00479 RooNLLVar b_nll("b_nll","b_nll",*fBModel,*dataToTest,RooFit::CloneData(false),RooFit::Extended());
00480 double m2lnQ = -2*(b_nll.getVal()-sb_nll.getVal());
00481 if ( hypoTested==0 ) sbVals.push_back(m2lnQ);
00482 else bVals.push_back(m2lnQ);
00483 } else {
00484 RooNLLVar sb_nll("sb_nll","sb_nll",*fSbModel,*dataToTest,RooFit::CloneData(false));
00485 RooNLLVar b_nll("b_nll","b_nll",*fBModel,*dataToTest,RooFit::CloneData(false));
00486 double m2lnQ = -2*(b_nll.getVal()-sb_nll.getVal());
00487 if ( hypoTested==0 ) sbVals.push_back(m2lnQ);
00488 else bVals.push_back(m2lnQ);
00489 }
00490 }
00491 }
00492
00493
00494 delete sbData;
00495 delete bData;
00496
00497
00498 if (fTestStatisticsIdx == 3) {
00499 RooArgSet * sbparams = fSbModel->getParameters(*fObservables);
00500 if (sbparams) {
00501 assert(originalSbParams.getSize() == sbparams->getSize());
00502 *sbparams = originalSbParams;
00503 delete sbparams;
00504 }
00505 RooArgSet * bparams = fBModel->getParameters(*fObservables);
00506 if (bparams) {
00507 assert(originalBParams.getSize() == bparams->getSize());
00508 *bparams = originalBParams;
00509 delete bparams;
00510 }
00511 }
00512
00513
00514
00515 }
00516
00517
00518
00519 if (usePriors && nParameters>0) {
00520 for (int iParameter=0; iParameter<nParameters; iParameter++) {
00521 RooRealVar* oneParam = (RooRealVar*) parametersList.at(iParameter);
00522 oneParam->setVal(parameterValues[iParameter]);
00523 }
00524 }
00525
00526 return;
00527 }
00528
00529
00530
00531 void HybridCalculatorOriginal::PrintMore(const char* options) const
00532 {
00533
00534
00535 if (fSbModel) {
00536 std::cout << "Signal plus background model:\n";
00537 fSbModel->Print(options);
00538 }
00539
00540 if (fBModel) {
00541 std::cout << "\nBackground model:\n";
00542 fBModel->Print(options);
00543 }
00544
00545 if (fObservables) {
00546 std::cout << "\nObservables:\n";
00547 fObservables->Print(options);
00548 }
00549
00550 if (fNuisanceParameters) {
00551 std::cout << "\nParameters being integrated:\n";
00552 fNuisanceParameters->Print(options);
00553 }
00554
00555 if (fPriorPdf) {
00556 std::cout << "\nPrior PDF model for integration:\n";
00557 fPriorPdf->Print(options);
00558 }
00559
00560 return;
00561 }
00562
00563
00564
00565 HybridResult* HybridCalculatorOriginal::GetHypoTest() const {
00566
00567
00568
00569 if (!DoCheckInputs()) return 0;
00570 RooAbsData * treeData = dynamic_cast<RooAbsData *> (fData);
00571 if (!treeData) {
00572 std::cerr << "Error in HybridCalculatorOriginal::GetHypoTest - invalid data type - return NULL" << std::endl;
00573 return 0;
00574 }
00575 bool usePrior = (fUsePriorPdf && fPriorPdf );
00576 return Calculate( *treeData, fNToys, usePrior);
00577 }
00578
00579
00580 bool HybridCalculatorOriginal::DoCheckInputs() const {
00581 if (!fData) {
00582 std::cerr << "Error in HybridCalculatorOriginal - data have not been set" << std::endl;
00583 return false;
00584 }
00585
00586
00587 if (!fObservables && fData->get() ) fObservables = new RooArgList( *fData->get() );
00588 if (!fObservables) {
00589 std::cerr << "Error in HybridCalculatorOriginal - no observables" << std::endl;
00590 return false;
00591 }
00592
00593 if (!fSbModel) {
00594 std::cerr << "Error in HybridCalculatorOriginal - S+B pdf has not been set " << std::endl;
00595 return false;
00596 }
00597
00598 if (!fBModel) {
00599 std::cerr << "Error in HybridCalculatorOriginal - B pdf has not been set" << std::endl;
00600 return false;
00601 }
00602 if (fUsePriorPdf && !fNuisanceParameters) {
00603 std::cerr << "Error in HybridCalculatorOriginal - nuisance parameters have not been set " << std::endl;
00604 return false;
00605 }
00606 if (fUsePriorPdf && !fPriorPdf) {
00607 std::cerr << "Error in HybridCalculatorOriginal - prior pdf has not been set " << std::endl;
00608 return false;
00609 }
00610 return true;
00611 }
00612
00613