//put very small data entries in a binned dataset to avoid unphysical pdfs, specifically for H->ZZ->4l RooDataSet* makeData(RooDataSet* orig, RooSimultaneous* simPdf, const RooArgSet* observables, RooRealVar* firstPOI, double mass, double& mu_min) { double max_soverb = 0; mu_min = -10e9; map<string, RooDataSet*> data_map; firstPOI->setVal(0); RooCategory* cat = (RooCategory*)&simPdf->indexCat(); TList* datalist = orig->split(*(RooAbsCategory*)cat, true); TIterator* dataItr = datalist->MakeIterator(); RooAbsData* ds; RooRealVar* weightVar = new RooRealVar("weightVar","weightVar",1); while ((ds = (RooAbsData*)dataItr->Next())) { string typeName(ds->GetName()); cat->setLabel(typeName.c_str()); RooAbsPdf* pdf = simPdf->getPdf(typeName.c_str()); cout << "pdf: " << pdf << endl; RooArgSet* obs = pdf->getObservables(observables); cout << "obs: " << obs << endl; RooArgSet obsAndWeight(*obs, *weightVar); obsAndWeight.add(*cat); stringstream datasetName; datasetName << "newData_" << typeName; RooDataSet* thisData = new RooDataSet(datasetName.str().c_str(),datasetName.str().c_str(), obsAndWeight, WeightVar(*weightVar)); RooRealVar* firstObs = (RooRealVar*)obs->first(); //int ibin = 0; int nrEntries = ds->numEntries(); for (int ib=0;ib<nrEntries;ib++) { const RooArgSet* event = ds->get(ib); const RooRealVar* thisObs = (RooRealVar*)event->find(firstObs->GetName()); firstObs->setVal(thisObs->getVal()); firstPOI->setVal(0); double b = pdf->expectedEvents(*firstObs)*pdf->getVal(obs); firstPOI->setVal(1); double s = pdf->expectedEvents(*firstObs)*pdf->getVal(obs) - b; if (s > 0) { mu_min = max(mu_min, -b/s); double soverb = s/b; if (soverb > max_soverb) { max_soverb = soverb; cout << "Found new max s/b: " << soverb << " in pdf " << pdf->GetName() << " at m = " << thisObs->getVal() << endl; } } if (b == 0 && s != 0) { cout << "Expecting non-zero signal and zero bg at m=" << firstObs->getVal() << " in pdf " << pdf->GetName() << endl; } if (s+b <= 0) { cout << "expecting zero" << endl; continue; } double weight = ds->weight(); if ((typeName.find("ATLAS_H_4mu") != string::npos || typeName.find("ATLAS_H_4e") != string::npos || typeName.find("ATLAS_H_2mu2e") != string::npos || typeName.find("ATLAS_H_2e2mu") != string::npos) && fabs(firstObs->getVal() - mass) < 10 && weight == 0) { cout << "adding event: " << firstObs->getVal() << endl; thisData->add(*event, pow(10., -9.)); } else { //weight = max(pow(10.0, -9), weight); thisData->add(*event, weight); } } data_map[string(ds->GetName())] = (RooDataSet*)thisData; } RooDataSet* newData = new RooDataSet("newData","newData",RooArgSet(*observables, *weightVar), Index(*cat), Import(data_map), WeightVar(*weightVar)); orig->Print(); newData->Print(); //newData->tree()->Scan("*"); return newData; }
void StandardHistFactoryPlotsWithCategories(const char* infile = "", const char* workspaceName = "combined", const char* modelConfigName = "ModelConfig", const char* dataName = "obsData"){ double nSigmaToVary=5.; double muVal=0; bool doFit=false; // ------------------------------------------------------- // First part is just to access a user-defined file // or create the standard example file if it doesn't exist const char* filename = ""; if (!strcmp(infile,"")) { filename = "results/example_combined_GaussExample_model.root"; bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code // if file does not exists generate with histfactory if (!fileExist) { #ifdef _WIN32 cout << "HistFactory file cannot be generated on Windows - exit" << endl; return; #endif // Normally this would be run on the command line cout <<"will run standard hist2workspace example"<<endl; gROOT->ProcessLine(".! prepareHistFactory ."); gROOT->ProcessLine(".! hist2workspace config/example.xml"); cout <<"\n\n---------------------"<<endl; cout <<"Done creating example input"<<endl; cout <<"---------------------\n\n"<<endl; } } else filename = infile; // Try to open the file TFile *file = TFile::Open(filename); // if input file was specified byt not found, quit if(!file ){ cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl; return; } // ------------------------------------------------------- // Tutorial starts here // ------------------------------------------------------- // get the workspace out of the file RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName); if(!w){ cout <<"workspace not found" << endl; return; } // get the modelConfig out of the file ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName); // get the modelConfig out of the file RooAbsData* data = w->data(dataName); // make sure ingredients are found if(!data || !mc){ w->Print(); cout << "data or ModelConfig was not found" <<endl; return; } // ------------------------------------------------------- // now use the profile inspector RooRealVar* obs = (RooRealVar*)mc->GetObservables()->first(); TList* list = new TList(); RooRealVar * firstPOI = dynamic_cast<RooRealVar*>(mc->GetParametersOfInterest()->first()); firstPOI->setVal(muVal); // firstPOI->setConstant(); if(doFit){ mc->GetPdf()->fitTo(*data); } // ------------------------------------------------------- mc->GetNuisanceParameters()->Print("v"); int nPlotsMax = 1000; cout <<" check expectedData by category"<<endl; RooDataSet* simData=NULL; RooSimultaneous* simPdf = NULL; if(strcmp(mc->GetPdf()->ClassName(),"RooSimultaneous")==0){ cout <<"Is a simultaneous PDF"<<endl; simPdf = (RooSimultaneous *)(mc->GetPdf()); } else { cout <<"Is not a simultaneous PDF"<<endl; } if(doFit) { RooCategory* channelCat = (RooCategory*) (&simPdf->indexCat()); TIterator* iter = channelCat->typeIterator() ; RooCatType* tt = NULL; tt=(RooCatType*) iter->Next(); RooAbsPdf* pdftmp = ((RooSimultaneous*)mc->GetPdf())->getPdf(tt->GetName()) ; RooArgSet* obstmp = pdftmp->getObservables(*mc->GetObservables()) ; obs = ((RooRealVar*)obstmp->first()); RooPlot* frame = obs->frame(); cout <<Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())<<endl; cout << tt->GetName() << " " << channelCat->getLabel() <<endl; data->plotOn(frame,MarkerSize(1),Cut(Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())),DataError(RooAbsData::None)); Double_t normCount = data->sumEntries(Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())) ; pdftmp->plotOn(frame,LineWidth(2.),Normalization(normCount,RooAbsReal::NumEvent)) ; frame->Draw(); cout <<"expected events = " << mc->GetPdf()->expectedEvents(*data->get()) <<endl; return; } int nPlots=0; if(!simPdf){ TIterator* it = mc->GetNuisanceParameters()->createIterator(); RooRealVar* var = NULL; while( (var = (RooRealVar*) it->Next()) != NULL){ RooPlot* frame = obs->frame(); frame->SetYTitle(var->GetName()); data->plotOn(frame,MarkerSize(1)); var->setVal(0); mc->GetPdf()->plotOn(frame,LineWidth(1.)); var->setVal(1); mc->GetPdf()->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(1)); var->setVal(-1); mc->GetPdf()->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(1)); list->Add(frame); var->setVal(0); } } else { RooCategory* channelCat = (RooCategory*) (&simPdf->indexCat()); // TIterator* iter = simPdf->indexCat().typeIterator() ; TIterator* iter = channelCat->typeIterator() ; RooCatType* tt = NULL; while(nPlots<nPlotsMax && (tt=(RooCatType*) iter->Next())) { cout << "on type " << tt->GetName() << " " << endl; // Get pdf associated with state from simpdf RooAbsPdf* pdftmp = simPdf->getPdf(tt->GetName()) ; // Generate observables defined by the pdf associated with this state RooArgSet* obstmp = pdftmp->getObservables(*mc->GetObservables()) ; // obstmp->Print(); obs = ((RooRealVar*)obstmp->first()); TIterator* it = mc->GetNuisanceParameters()->createIterator(); RooRealVar* var = NULL; while(nPlots<nPlotsMax && (var = (RooRealVar*) it->Next())){ TCanvas* c2 = new TCanvas("c2"); RooPlot* frame = obs->frame(); frame->SetName(Form("frame%d",nPlots)); frame->SetYTitle(var->GetName()); cout <<Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())<<endl; cout << tt->GetName() << " " << channelCat->getLabel() <<endl; data->plotOn(frame,MarkerSize(1),Cut(Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())),DataError(RooAbsData::None)); Double_t normCount = data->sumEntries(Form("%s==%s::%s",channelCat->GetName(),channelCat->GetName(),tt->GetName())) ; if(strcmp(var->GetName(),"Lumi")==0){ cout <<"working on lumi"<<endl; var->setVal(w->var("nominalLumi")->getVal()); var->Print(); } else{ var->setVal(0); } // w->allVars().Print("v"); // mc->GetNuisanceParameters()->Print("v"); // pdftmp->plotOn(frame,LineWidth(2.)); // mc->GetPdf()->plotOn(frame,LineWidth(2.),Slice(*channelCat,tt->GetName()),ProjWData(*data)); //pdftmp->plotOn(frame,LineWidth(2.),Slice(*channelCat,tt->GetName()),ProjWData(*data)); normCount = pdftmp->expectedEvents(*obs); pdftmp->plotOn(frame,LineWidth(2.),Normalization(normCount,RooAbsReal::NumEvent)) ; if(strcmp(var->GetName(),"Lumi")==0){ cout <<"working on lumi"<<endl; var->setVal(w->var("nominalLumi")->getVal()+0.05); var->Print(); } else{ var->setVal(nSigmaToVary); } // pdftmp->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2)); // mc->GetPdf()->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2.),Slice(*channelCat,tt->GetName()),ProjWData(*data)); //pdftmp->plotOn(frame,LineColor(kRed),LineStyle(kDashed),LineWidth(2.),Slice(*channelCat,tt->GetName()),ProjWData(*data)); normCount = pdftmp->expectedEvents(*obs); pdftmp->plotOn(frame,LineWidth(2.),LineColor(kRed),LineStyle(kDashed),Normalization(normCount,RooAbsReal::NumEvent)) ; if(strcmp(var->GetName(),"Lumi")==0){ cout <<"working on lumi"<<endl; var->setVal(w->var("nominalLumi")->getVal()-0.05); var->Print(); } else{ var->setVal(-nSigmaToVary); } // pdftmp->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2)); // mc->GetPdf()->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2),Slice(*channelCat,tt->GetName()),ProjWData(*data)); //pdftmp->plotOn(frame,LineColor(kGreen),LineStyle(kDashed),LineWidth(2),Slice(*channelCat,tt->GetName()),ProjWData(*data)); normCount = pdftmp->expectedEvents(*obs); pdftmp->plotOn(frame,LineWidth(2.),LineColor(kGreen),LineStyle(kDashed),Normalization(normCount,RooAbsReal::NumEvent)) ; // set them back to normal if(strcmp(var->GetName(),"Lumi")==0){ cout <<"working on lumi"<<endl; var->setVal(w->var("nominalLumi")->getVal()); var->Print(); } else{ var->setVal(0); } list->Add(frame); // quit making plots ++nPlots; frame->Draw(); c2->SaveAs(Form("%s_%s_%s.pdf",tt->GetName(),obs->GetName(),var->GetName())); delete c2; } } } // ------------------------------------------------------- // now make plots TCanvas* c1 = new TCanvas("c1","ProfileInspectorDemo",800,200); if(list->GetSize()>4){ double n = list->GetSize(); int nx = (int)sqrt(n) ; int ny = TMath::CeilNint(n/nx); nx = TMath::CeilNint( sqrt(n) ); c1->Divide(ny,nx); } else c1->Divide(list->GetSize()); for(int i=0; i<list->GetSize(); ++i){ c1->cd(i+1); list->At(i)->Draw(); } }