CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/HLTriggerOffline/HeavyFlavor/src/HeavyFlavorHarvesting.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/Frameworkfwd.h"
00002 #include "FWCore/Framework/interface/EDAnalyzer.h"
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/MakerMacros.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "FWCore/ServiceRegistry/interface/Service.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 #include "FWCore/Utilities/interface/InputTag.h"
00009 #include "DQMServices/Core/interface/DQMStore.h"
00010 #include "DQMServices/Core/interface/MonitorElement.h"
00011 
00012 #include "TH1F.h"
00013 #include "TH2F.h"
00014 #include "TH3F.h"
00015 #include "RVersion.h"
00016 
00017 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,27,0)
00018 #include "TEfficiency.h"
00019 #else
00020 #include "TGraphAsymmErrors.h"
00021 #endif
00022 using namespace edm;
00023 using namespace std;
00024 
00025 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,27,0)
00026 class HeavyFlavorHarvesting : public edm::EDAnalyzer { //, public TGraphAsymmErrors{
00027 #else
00028 class HeavyFlavorHarvesting : public edm::EDAnalyzer , public TGraphAsymmErrors{
00029 #endif
00030   public:
00031     HeavyFlavorHarvesting(const edm::ParameterSet& pset);
00032     virtual ~HeavyFlavorHarvesting();
00033     virtual void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) {};
00034     virtual void endRun(const edm::Run &, const edm::EventSetup &);
00035   private:
00036     void calculateEfficiency(const ParameterSet& pset);
00037     void calculateEfficiency1D( TH1* num, TH1* den, string name );
00038     void calculateEfficiency2D( TH2F* num, TH2F* den, string name );
00039     DQMStore * dqmStore;
00040     string myDQMrootFolder;
00041     const VParameterSet efficiencies;
00042 };
00043 
00044 HeavyFlavorHarvesting::HeavyFlavorHarvesting(const edm::ParameterSet& pset):
00045   myDQMrootFolder( pset.getUntrackedParameter<string>("MyDQMrootFolder") ),
00046   efficiencies( pset.getUntrackedParameter<VParameterSet>("Efficiencies") )
00047 {
00048 }
00049 
00050 void HeavyFlavorHarvesting::endRun(const edm::Run &, const edm::EventSetup &){
00051   dqmStore = Service<DQMStore>().operator->();
00052   if( !dqmStore ){
00053     LogError("HLTriggerOfflineHeavyFlavor") << "Could not find DQMStore service\n";
00054     return;
00055   }
00056   for(VParameterSet::const_iterator pset = efficiencies.begin(); pset!=efficiencies.end(); pset++){
00057     calculateEfficiency(*pset);
00058   }
00059 }
00060   
00061 void HeavyFlavorHarvesting::calculateEfficiency(const ParameterSet& pset){
00062 //get hold of numerator and denominator histograms
00063   vector<string> numDenEffMEnames = pset.getUntrackedParameter<vector<string> >("NumDenEffMEnames");
00064   if(numDenEffMEnames.size()!=3){
00065     LogDebug("HLTriggerOfflineHeavyFlavor") << "NumDenEffMEnames must have three names"<<endl;
00066     return;
00067   }
00068   string denMEname = myDQMrootFolder+"/"+numDenEffMEnames[1];
00069   string numMEname = myDQMrootFolder+"/"+numDenEffMEnames[0];
00070   MonitorElement *denME = dqmStore->get(denMEname);
00071   MonitorElement *numME = dqmStore->get(numMEname);
00072   if(denME==0 || numME==0){
00073     LogDebug("HLTriggerOfflineHeavyFlavor") << "Could not find MEs: "<<denMEname<<" or "<<numMEname<<endl;
00074     return;
00075   }
00076   TH1 *den = denME->getTH1();
00077   TH1 *num = numME->getTH1();
00078   //check compatibility of the histograms  
00079   if( den->GetNbinsX() != num->GetNbinsX() || den->GetNbinsY() != num->GetNbinsY() ||  den->GetNbinsZ() != num->GetNbinsZ() ){
00080     LogDebug("HLTriggerOfflineHeavyFlavor") << "Monitoring elements "<<numMEname<<" and "<<denMEname<<"are incompatible"<<endl;
00081     return;
00082   }
00083   //figure out the directory and efficiency name  
00084   string effName = numDenEffMEnames[2];
00085   string effDir = myDQMrootFolder;
00086   string::size_type slashPos = effName.rfind('/');
00087   if ( string::npos != slashPos ) {
00088     effDir += "/"+effName.substr(0, slashPos);
00089     effName.erase(0, slashPos+1);
00090   }
00091   dqmStore->setCurrentFolder(effDir);
00092   //calculate the efficiencies
00093   int dimensions = num->GetDimension();
00094   if(dimensions==1){
00095     calculateEfficiency1D( num, den, effName );
00096   }else if(dimensions==2){
00097     calculateEfficiency2D( (TH2F*)num, (TH2F*)den, effName );
00098     TH1D* numX = ((TH2F*)num)->ProjectionX(); 
00099     TH1D* denX = ((TH2F*)den)->ProjectionX();
00100     calculateEfficiency1D( numX, denX, effName+"X" );
00101     delete numX;
00102     delete denX;
00103     TH1D* numY = ((TH2F*)num)->ProjectionY();
00104     TH1D* denY = ((TH2F*)den)->ProjectionY();
00105     calculateEfficiency1D( numY, denY, effName+"Y" );
00106     delete numY;
00107     delete denY;
00108   }else{
00109     return;
00110   }
00111 }
00112 
00113 void HeavyFlavorHarvesting::calculateEfficiency1D( TH1* num, TH1* den, string effName ){
00114   TProfile* eff;
00115   if(num->GetXaxis()->GetXbins()->GetSize()==0){
00116     eff = new TProfile(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXmin(),num->GetXaxis()->GetXmax());
00117   }else{
00118     eff = new TProfile(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXbins()->GetArray());
00119   }
00120   eff->SetTitle(effName.c_str());
00121   eff->SetXTitle( num->GetXaxis()->GetTitle() );
00122   eff->SetYTitle("Efficiency");
00123   eff->SetOption("PE");
00124   eff->SetLineColor(2);
00125   eff->SetLineWidth(2);
00126   eff->SetMarkerStyle(20);
00127   eff->SetMarkerSize(0.8);
00128   eff->GetYaxis()->SetRangeUser(-0.001,1.001);
00129   eff->SetStats(kFALSE);
00130   for(int i=1;i<=num->GetNbinsX();i++){
00131     double e, low, high;
00132 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,27,0)
00133     if (int(den->GetBinContent(i))>0.) e= double(num->GetBinContent(i))/double(den->GetBinContent(i));
00134     else e=0.;
00135     low=TEfficiency::Wilson((double)den->GetBinContent(i),(double)num->GetBinContent(i),0.683,false);
00136     high=TEfficiency::Wilson((double)den->GetBinContent(i),(double)num->GetBinContent(i),0.683,true);
00137 #else
00138     Efficiency( (double)num->GetBinContent(i), (double)den->GetBinContent(i), 0.683, e, low, high );
00139 #endif
00140     double err = e-low>high-e ? e-low : high-e;
00141     //here is the trick to store info in TProfile:
00142     eff->SetBinContent( i, e );
00143     eff->SetBinEntries( i, 1 );
00144     eff->SetBinError( i, sqrt(e*e+err*err) );
00145   }
00146   dqmStore->bookProfile(effName,eff);
00147   delete eff;
00148 }
00149 
00150 void HeavyFlavorHarvesting::calculateEfficiency2D( TH2F* num, TH2F* den, string effName ){
00151   TProfile2D* eff;
00152   if(num->GetXaxis()->GetXbins()->GetSize()==0 && num->GetYaxis()->GetXbins()->GetSize()==0){
00153     eff = new TProfile2D(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXmin(),num->GetXaxis()->GetXmax(),num->GetYaxis()->GetNbins(),num->GetYaxis()->GetXmin(),num->GetYaxis()->GetXmax());
00154   }else if(num->GetXaxis()->GetXbins()->GetSize()!=0 && num->GetYaxis()->GetXbins()->GetSize()==0){
00155     eff = new TProfile2D(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXbins()->GetArray(),num->GetYaxis()->GetNbins(),num->GetYaxis()->GetXmin(),num->GetYaxis()->GetXmax());
00156   }else if(num->GetXaxis()->GetXbins()->GetSize()==0 && num->GetYaxis()->GetXbins()->GetSize()!=0){
00157     eff = new TProfile2D(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXmin(),num->GetXaxis()->GetXmax(),num->GetYaxis()->GetNbins(),num->GetYaxis()->GetXbins()->GetArray());
00158   }else{
00159     eff = new TProfile2D(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXbins()->GetArray(),num->GetYaxis()->GetNbins(),num->GetYaxis()->GetXbins()->GetArray());
00160   }
00161   eff->SetTitle(effName.c_str());
00162   eff->SetXTitle( num->GetXaxis()->GetTitle() );
00163   eff->SetYTitle( num->GetYaxis()->GetTitle() );
00164   eff->SetZTitle("Efficiency");
00165   eff->SetOption("colztexte");
00166   eff->GetZaxis()->SetRangeUser(-0.001,1.001);
00167   eff->SetStats(kFALSE);
00168   for(int i=0;i<num->GetSize();i++){
00169     double e, low, high;
00170 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,27,0)
00171     if (int(den->GetBinContent(i))>0.) e= double(num->GetBinContent(i))/double(den->GetBinContent(i));
00172     else e=0.;
00173     low=TEfficiency::Wilson((double)den->GetBinContent(i),(double)num->GetBinContent(i),0.683,false);
00174     high=TEfficiency::Wilson((double)den->GetBinContent(i),(double)num->GetBinContent(i),0.683,true);
00175 #else
00176     Efficiency( (double)num->GetBinContent(i), (double)den->GetBinContent(i), 0.683, e, low, high );
00177 #endif
00178     double err = e-low>high-e ? e-low : high-e;
00179     //here is the trick to store info in TProfile:
00180     eff->SetBinContent( i, e );
00181     eff->SetBinEntries( i, 1 );
00182     eff->SetBinError( i, sqrt(e*e+err*err) );
00183   }
00184   dqmStore->bookProfile2D(effName,eff);
00185   delete eff;
00186 }
00187 
00188 HeavyFlavorHarvesting::~HeavyFlavorHarvesting(){
00189 }
00190 
00191 //define this as a plug-in
00192 DEFINE_FWK_MODULE(HeavyFlavorHarvesting);