CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/DQM/DataScouting/plugins/DiJetVarAnalyzer.cc

Go to the documentation of this file.
00001 #include "DQM/DataScouting/plugins/DiJetVarAnalyzer.h"
00002 
00003 #include "DataFormats/JetReco/interface/CaloJet.h"
00004 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00005 
00006 #include "DataFormats/METReco/interface/CaloMET.h"
00007 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00008 
00009 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00010 #include "DataFormats/EgammaCandidates/interface/ElectronFwd.h"
00011 
00012 #include "DataFormats/MuonReco/interface/Muon.h"
00013 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00014 
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 
00017 #include "HLTrigger/HLTcore/interface/TriggerExpressionData.h"
00018 #include "HLTrigger/HLTcore/interface/TriggerExpressionEvaluator.h"
00019 #include "HLTrigger/HLTcore/interface/TriggerExpressionParser.h"
00020 
00021 #include <cmath>
00022 
00023 //------------------------------------------------------------------------------
00024 // A simple constructor which takes as inoput only the name of the PF jet collection
00025 DiJetVarAnalyzer::DiJetVarAnalyzer( const edm::ParameterSet & conf ):
00026   ScoutingAnalyzerBase(conf),
00027   jetCollectionTag_        (conf.getUntrackedParameter<edm::InputTag>("jetCollectionTag")),
00028   //dijetVarCollectionTag_   (conf.getUntrackedParameter<edm::InputTag>("dijetVarCollectionTag")),
00029   widejetsCollectionTag_   (conf.getUntrackedParameter<edm::InputTag>("widejetsCollectionTag")),
00030   metCollectionTag_        (conf.getUntrackedParameter<edm::InputTag>("metCollectionTag")),
00031   metCleanCollectionTag_   (conf.getUntrackedParameter<edm::InputTag>("metCleanCollectionTag")),
00032   numwidejets_             (conf.getParameter<unsigned int>("numwidejets")),
00033   etawidejets_             (conf.getParameter<double>("etawidejets")),
00034   ptwidejets_              (conf.getParameter<double>("ptwidejets")),
00035   detawidejets_            (conf.getParameter<double>("detawidejets")),
00036   dphiwidejets_            (conf.getParameter<double>("dphiwidejets")),
00037   maxEMfraction_           (conf.getParameter<double>("maxEMfraction")),
00038   maxHADfraction_          (conf.getParameter<double>("maxHADfraction")),
00039   HLTpathMain_             (triggerExpression::parse( conf.getParameter<std::string>("HLTpathMain") )),
00040   HLTpathMonitor_          (triggerExpression::parse( conf.getParameter<std::string>("HLTpathMonitor") )),
00041   triggerConfiguration_           (conf.getParameterSet("triggerConfiguration"))
00042 {
00043 }
00044 
00045 //------------------------------------------------------------------------------
00046 // Nothing to destroy: the DQM service thinks about everything
00047 DiJetVarAnalyzer::~DiJetVarAnalyzer(){
00048   delete HLTpathMain_;
00049   delete HLTpathMonitor_;
00050 }
00051 
00052 //------------------------------------------------------------------------------
00053 void DiJetVarAnalyzer::beginRun( const edm::Run & iRun, const edm::EventSetup & iSetup){
00054 }
00055 
00056 //------------------------------------------------------------------------------
00057 // Usual analyze method
00058 void DiJetVarAnalyzer::analyze( const edm::Event & iEvent, const edm::EventSetup & c ){
00059   
00060   using namespace std;
00061   using namespace edm;
00062   using namespace reco;
00063   
00064   // ## Get jet collection
00065   edm::Handle<reco::CaloJetCollection> calojets_handle;
00066   iEvent.getByLabel(jetCollectionTag_,calojets_handle);
00067 
00068   // Loop over all the selected jets ( defined at DQM/DataScouting/python/dijetScouting_cff.py )  
00069   double thisHT = 0.;
00070   for(reco::CaloJetCollection::const_iterator it = calojets_handle->begin(); it != calojets_handle->end(); ++it)
00071     {
00072       //cout << "== jet: " << it->pt() << " " << it->eta() << " " << it->phi() << endl;
00073       m_selJets_pt->Fill( it->pt() );
00074       m_selJets_eta->Fill( it->eta() );
00075       m_selJets_phi->Fill( it->phi() );
00076       m_selJets_hadEnergyFraction->Fill( it->energyFractionHadronic() );
00077       m_selJets_emEnergyFraction->Fill( it->emEnergyFraction() );
00078       m_selJets_towersArea->Fill( it->towersArea() );
00079 
00080       thisHT += it->pt();
00081     }
00082 
00083   // HT
00084   m_HT_inclusive->Fill(thisHT);      
00085   
00086   // ## Get widejets 
00087   edm::Handle< vector<math::PtEtaPhiMLorentzVector> > widejets_handle;
00088   iEvent.getByLabel (widejetsCollectionTag_,widejets_handle);
00089   
00090   TLorentzVector wj1;
00091   TLorentzVector wj2;
00092   TLorentzVector wdijet;
00093 
00094   double MJJWide = -1;
00095   double DeltaEtaJJWide = -1;
00096   double DeltaPhiJJWide = -1;
00097 
00098   if( widejets_handle->size() >= 2 )
00099     {
00100       wj1.SetPtEtaPhiM(widejets_handle->at(0).pt(),
00101                        widejets_handle->at(0).eta(),
00102                        widejets_handle->at(0).phi(),
00103                        widejets_handle->at(0).mass()
00104                        );
00105       wj2.SetPtEtaPhiM(widejets_handle->at(1).pt(),
00106                        widejets_handle->at(1).eta(),
00107                        widejets_handle->at(1).phi(),
00108                        widejets_handle->at(1).mass()
00109                        );
00110 
00111       wdijet = wj1 + wj2;
00112 
00113       MJJWide = wdijet.M();
00114       DeltaEtaJJWide = fabs(wj1.Eta()-wj2.Eta());
00115       DeltaPhiJJWide = fabs(wj1.DeltaPhi(wj2));
00116 
00117       //       cout << "== j1 wide: " << wj1.Pt() << " " << wj1.Eta() << " " << wj1.Phi() << " " << wj1.M() << endl;
00118       //       cout << "== j2 wide: " << wj2.Pt() << " " << wj2.Eta() << " " << wj2.Phi() << " " << wj2.M() << endl;
00119       //       cout << "== MJJWide: " << MJJWide << endl;
00120       //       cout << "== DeltaEtaJJWide: " << DeltaEtaJJWide << endl;
00121       //       cout << "== DeltaPhiJJWide: " << DeltaPhiJJWide << endl;
00122     }
00123 
00124   // ## Get met collection
00125 
00126   // met
00127   edm::Handle<reco::CaloMETCollection> calomet_handle;
00128   iEvent.getByLabel(metCollectionTag_,calomet_handle);
00129 
00130   // met cleaned
00131   edm::Handle<reco::CaloMETCollection> calometClean_handle;
00132   iEvent.getByLabel(metCleanCollectionTag_,calometClean_handle);
00133  
00134   if( calomet_handle.isValid() && calometClean_handle.isValid() )
00135     {
00136       //       std::cout << "---" << std::endl;
00137       //       std::cout << "== calomet: " << (calomet_handle->front()).pt() << " " << (calomet_handle->front()).phi() << std::endl;
00138       //       std::cout << "== calometClean: " << (calometClean_handle->front()).pt() << " " << (calometClean_handle->front()).phi() << std::endl;
00139       //       std::cout << "== calomet - calometClean: " << (calomet_handle->front()).pt() - (calometClean_handle->front()).pt() << std::endl;
00140       //       std::cout << "---" << std::endl;
00141       m_metCases->Fill(0);
00142       m_metDiff->Fill( (calomet_handle->front()).pt() - (calometClean_handle->front()).pt() );
00143       m_metVSmetclean->Fill( (calomet_handle->front()).pt() , (calometClean_handle->front()).pt() );      
00144     }
00145   else if( calomet_handle.isValid() && !calometClean_handle.isValid() )
00146     {
00147       m_metCases->Fill(1);
00148       m_metCaseNoMetClean->Fill((calomet_handle->front()).pt());
00149     }
00150   else if( !calomet_handle.isValid() && calometClean_handle.isValid() )
00151     {
00152       m_metCases->Fill(2);
00153     }
00154 
00155   // ## Event Selection
00156   bool pass_nocut=false;
00157   bool pass_twowidejets=false;
00158   bool pass_etaptcuts=false;
00159   bool pass_deta=false;
00160   bool pass_JetIDtwojets=true;
00161   bool pass_dphi=false;
00162   bool pass_metFilter=true;
00163   //--
00164   bool pass_deta_L4=false;
00165   bool pass_deta_L3=false;
00166   bool pass_deta_L2=false;
00167 
00168   bool pass_fullsel_NOdeta=false;
00169   bool pass_fullsel_detaL4=false;
00170   bool pass_fullsel_detaL3=false;
00171   bool pass_fullsel_detaL2=false;
00172   bool pass_fullsel=false;
00173 
00174   // No cut
00175   pass_nocut=true;
00176 
00177   // Two wide jets
00178   if( widejets_handle->size() >= numwidejets_ )
00179     {
00180       // Two wide jets
00181       pass_twowidejets=true;
00182 
00183       // Eta/pt cuts
00184       if( fabs(wj1.Eta())<etawidejets_ && wj1.Pt()>ptwidejets_ 
00185           &&
00186           fabs(wj2.Eta())<etawidejets_ && wj2.Pt()>ptwidejets_
00187           )
00188         {
00189           pass_etaptcuts=true;
00190         }
00191       
00192       // Deta cut
00193       if( DeltaEtaJJWide < detawidejets_ )
00194         pass_deta=true;
00195 
00196       // Dphi cut
00197       if( DeltaPhiJJWide > dphiwidejets_ )
00198         pass_dphi=true;
00199 
00200       // Other Deta cuts
00201       if( DeltaEtaJJWide < 4.0 )
00202         pass_deta_L4=true;
00203 
00204       if( DeltaEtaJJWide < 3.0 )
00205         pass_deta_L3=true;
00206 
00207       if( DeltaEtaJJWide < 2.0 )
00208         pass_deta_L2=true;
00209       
00210     }
00211 
00212   // Jet id two leading jets
00213 
00214   if( calojets_handle->size() >= numwidejets_ )
00215     {
00216       //   first jet
00217       reco::CaloJetCollection::const_iterator thisJet = calojets_handle->begin();
00218       //cout << "== thisJet1: " << thisJet->pt() << " " << thisJet->eta() << " " << thisJet->phi() << endl;
00219       if( thisJet->energyFractionHadronic()>maxHADfraction_ || thisJet->emEnergyFraction()>maxEMfraction_ )
00220         pass_JetIDtwojets=false;
00221 
00222       //   second jet
00223       ++thisJet;
00224       //cout << "== thisJet2: " << thisJet->pt() << " " << thisJet->eta() << " " << thisJet->phi() << endl;
00225       if( thisJet->energyFractionHadronic()>maxHADfraction_ || thisJet->emEnergyFraction()>maxEMfraction_ )
00226         pass_JetIDtwojets=false;
00227     }      
00228 
00229   // Met filter
00230   if( calomet_handle.isValid() && calometClean_handle.isValid() )
00231     {
00232       if( fabs ( (calomet_handle->front()).pt() - (calometClean_handle->front()).pt() ) > 0.1 )
00233         pass_metFilter=false;
00234     }
00235   
00236   // Full selection (no deta cut)
00237   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter )
00238     pass_fullsel_NOdeta=true;
00239 
00240   // Full selection (various deta cuts)
00241   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter && pass_deta_L4 )
00242     pass_fullsel_detaL4=true;
00243   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter && pass_deta_L3 )
00244     pass_fullsel_detaL3=true;
00245   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter && pass_deta_L2 )
00246     pass_fullsel_detaL2=true;
00247 
00248   // Full selection (default deta cut)
00249   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets && pass_dphi && pass_metFilter )
00250     pass_fullsel=true;
00251   
00252   // ## Fill Histograms 
00253 
00254   // Cut-flow plot
00255   if( pass_nocut )
00256     m_cutFlow->Fill(0);
00257   if( pass_nocut && pass_twowidejets )
00258     m_cutFlow->Fill(1);
00259   if( pass_nocut && pass_twowidejets && pass_etaptcuts )
00260     m_cutFlow->Fill(2);
00261   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta )
00262     m_cutFlow->Fill(3);
00263   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets )
00264     m_cutFlow->Fill(4);
00265   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets && pass_dphi )
00266     m_cutFlow->Fill(5);
00267   if( pass_fullsel )
00268     m_cutFlow->Fill(6);
00269 
00270   // After full selection
00271   if( pass_fullsel ) 
00272     {
00273       // 1D histograms
00274       m_MjjWide_finalSel->Fill(MJJWide);
00275       m_MjjWide_finalSel_varbin->Fill(MJJWide);
00276       m_DetajjWide_finalSel->Fill(DeltaEtaJJWide);
00277       m_DphijjWide_finalSel->Fill(DeltaPhiJJWide);      
00278       m_HT_finalSel->Fill(thisHT);      
00279     }      
00280 
00281   // After full selection (without "noise" filters)
00282   if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta )
00283     {
00284       m_MjjWide_finalSel_WithoutNoiseFilter->Fill(MJJWide);
00285       m_MjjWide_finalSel_WithoutNoiseFilter_varbin->Fill(MJJWide);
00286     }
00287 
00288   // After full selection (except DeltaEta cut)
00289   if( pass_fullsel_NOdeta )
00290     {
00291       // 1D histograms
00292       m_DetajjWide->Fill(DeltaEtaJJWide);
00293 
00294       if( DeltaEtaJJWide >= 0.0 && DeltaEtaJJWide < 0.5 )
00295         m_MjjWide_deta_0p0_0p5->Fill(MJJWide);
00296       if( DeltaEtaJJWide >= 0.5 && DeltaEtaJJWide < 1.0 )
00297         m_MjjWide_deta_0p5_1p0->Fill(MJJWide);
00298       if( DeltaEtaJJWide >= 1.0 && DeltaEtaJJWide < 1.5 )
00299         m_MjjWide_deta_1p0_1p5->Fill(MJJWide);
00300       if( DeltaEtaJJWide >= 1.5 && DeltaEtaJJWide < 2.0 )
00301         m_MjjWide_deta_1p5_2p0->Fill(MJJWide);
00302       if( DeltaEtaJJWide >= 2.0 && DeltaEtaJJWide < 2.5 )
00303         m_MjjWide_deta_2p0_2p5->Fill(MJJWide);
00304       if( DeltaEtaJJWide >= 2.5 && DeltaEtaJJWide < 3.0 )
00305         m_MjjWide_deta_2p5_3p0->Fill(MJJWide);
00306       if( DeltaEtaJJWide >= 3.0 )
00307         m_MjjWide_deta_3p0_inf->Fill(MJJWide);
00308       
00309       // 2D histograms
00310       m_DetajjVsMjjWide->Fill(MJJWide,DeltaEtaJJWide);            
00311       m_DetajjVsMjjWide_rebin->Fill(MJJWide,DeltaEtaJJWide);
00312     }
00313   
00314 
00315   // ## Get Trigger Info
00316 
00317   // HLT paths for DataScouting
00318   //  DST_HT250_v1
00319   //  DST_L1HTT_Or_L1MultiJet_v1
00320   //  DST_Mu5_HT250_v1
00321   //  DST_Ele8_CaloIdL_CaloIsoVL_TrkIdVL_TrkIsoVL_HT250_v1
00322 
00323   int HLTpathMain_fired    = -1;
00324   int HLTpathMonitor_fired = -1;
00325   
00326 
00327   if (HLTpathMain_ and HLTpathMonitor_ and triggerConfiguration_.setEvent(iEvent, c)) {
00328     // invalid HLT configuration, skip the processing
00329     
00330     // if the L1 or HLT configurations have changed, (re)initialize the filters (including during the first event)
00331     if (triggerConfiguration_.configurationUpdated()) {
00332       HLTpathMain_->init(triggerConfiguration_);
00333       HLTpathMonitor_->init(triggerConfiguration_);
00334       
00335       // log the expanded configuration
00336       // std::cout << "HLT selector configurations updated" << std::endl;
00337       // std::cout << "HLTpathMain:    " << *HLTpathMain_    << std::endl;
00338       // std::cout << "HLTpathMonitor: " << *HLTpathMonitor_ << std::endl;
00339     }
00340     
00341     HLTpathMain_fired    = (*HLTpathMain_)(triggerConfiguration_);
00342     HLTpathMonitor_fired = (*HLTpathMonitor_)(triggerConfiguration_);
00343     
00344     // The OR of the two should always be "1"
00345     // std::cout << *HLTpathMain_ << ": " << HLTpathMain_fired << " -- " << *HLTpathMonitor_ << ": " << HLTpathMonitor_fired << std::endl;
00346   }
00347   
00348   // ## Trigger Efficiency Curves
00349 
00350   //denominator - full sel NO deta cut
00351   if( pass_fullsel_NOdeta && HLTpathMonitor_fired == 1 )
00352     {
00353       m_MjjWide_den_NOdeta->Fill(MJJWide);
00354 
00355       //numerator  
00356       if( HLTpathMain_fired == 1)
00357         {
00358           m_MjjWide_num_NOdeta->Fill(MJJWide);
00359         }
00360     }
00361 
00362   //denominator - full sel deta < 4.0
00363   if( pass_fullsel_detaL4 && HLTpathMonitor_fired == 1 )
00364     {
00365       m_MjjWide_den_detaL4->Fill(MJJWide);
00366 
00367       //numerator  
00368       if( HLTpathMain_fired == 1)
00369         {
00370           m_MjjWide_num_detaL4->Fill(MJJWide);
00371         }
00372     }
00373 
00374   //denominator - full sel deta < 3.0
00375   if( pass_fullsel_detaL3 && HLTpathMonitor_fired == 1 )
00376     {
00377       m_MjjWide_den_detaL3->Fill(MJJWide);
00378 
00379       //numerator  
00380       if( HLTpathMain_fired == 1)
00381         {
00382           m_MjjWide_num_detaL3->Fill(MJJWide);
00383         }
00384     }
00385 
00386   //denominator - full sel deta < 2.0
00387   if( pass_fullsel_detaL2 && HLTpathMonitor_fired == 1 )
00388     {
00389       m_MjjWide_den_detaL2->Fill(MJJWide);
00390 
00391       //numerator  
00392       if( HLTpathMain_fired == 1)
00393         {
00394           m_MjjWide_num_detaL2->Fill(MJJWide);
00395         }
00396     }
00397   
00398   //denominator - full sel default deta cut (typically 1.3)
00399   if( pass_fullsel && HLTpathMonitor_fired == 1 )
00400     {
00401       m_MjjWide_den->Fill(MJJWide);
00402 
00403       //numerator  
00404       if( HLTpathMain_fired == 1)
00405         {
00406           m_MjjWide_num->Fill(MJJWide);
00407         }
00408     }
00409 
00410 
00411 }
00412 
00413 void DiJetVarAnalyzer::endRun( edm::Run const &, edm::EventSetup const & ){
00414 }
00415 
00416 //------------------------------------------------------------------------------
00417 // Function to book the Monitoring Elements.
00418 void DiJetVarAnalyzer::bookMEs(){
00419   
00420   // ==> TO BE UPDATED FOR sqrt(s)=8 TeV
00421   const int N_mass_bins=83;
00422   float massBins[N_mass_bins+1] = {1, 3, 6, 10, 16, 23, 31, 40, 50, 61, 74, 88, 103, 119, 137, 156, 176, 197, 220, 244, 270, 296, 325, 354, 386, 419, 453, 489, 526, 565, 606, 649, 693, 740, 788, 838, 890, 944, 1000, 1058, 1118, 1181, 1246, 1313, 1383, 1455, 1530, 1607, 1687, 1770, 1856, 1945, 2037, 2132, 2231, 2332, 2438, 2546, 2659, 2775, 2895, 3019, 3147, 3279, 3416, 3558, 3704, 3854, 4010, 4171, 4337, 4509, 4686, 4869, 5058, 5253, 5455, 5663, 5877, 6099, 6328, 6564, 6808, 7000};
00423 
00424   // 1D histograms
00425   m_cutFlow = bookH1withSumw2( "h1_cutFlow",
00426                                "Cut Flow",
00427                                7,0.,7.,
00428                                "Cut",
00429                                "Number of events"
00430                                );
00431   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(1,"No cut");
00432   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(2,"N(WideJets)>=2");
00433   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(3,"|#eta|<2.5 , pT>30");
00434   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(4,"|#Delta#eta|<1.3");
00435   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(5,"JetID");
00436   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(6,"|#Delta#phi|>#pi/3");
00437   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(7,"|met-metClean|>0.1");
00438 
00439   m_MjjWide_finalSel = bookH1withSumw2( "h1_MjjWide_finalSel",
00440                                         "M_{jj} WideJets (final selection)",
00441                                         8000,0.,8000.,
00442                                         "M_{jj} WideJets [GeV]",
00443                                         "Number of events"
00444                                         );
00445 
00446   m_MjjWide_finalSel_varbin = bookH1withSumw2BinArray( "h1_MjjWide_finalSel_varbin",
00447                                                        "M_{jj} WideJets (final selection)",
00448                                                        N_mass_bins, massBins,
00449                                                        "M_{jj} WideJets [GeV]",
00450                                                        "Number of events"
00451                                                        );
00452 
00453   m_MjjWide_finalSel_WithoutNoiseFilter = bookH1withSumw2( "h1_MjjWide_finalSel_WithoutNoiseFilter",
00454                                                            "M_{jj} WideJets (final selection, without noise filters)",
00455                                                            8000,0.,8000.,
00456                                                            "M_{jj} WideJets [GeV]",
00457                                                            "Number of events"
00458                                                            );
00459 
00460   m_MjjWide_finalSel_WithoutNoiseFilter_varbin = bookH1withSumw2BinArray( "h1_MjjWide_finalSel_WithoutNoiseFilter_varbin",
00461                                                                           "M_{jj} WideJets (final selection, without noise filters)",
00462                                                                           N_mass_bins, massBins,
00463                                                                           "M_{jj} WideJets [GeV]",
00464                                                                           "Number of events"
00465                                                                           );
00466   
00467   m_MjjWide_deta_0p0_0p5 = bookH1withSumw2( "h1_MjjWide_deta_0p0_0p5",
00468                                             "M_{jj} WideJets (0.0<=#Delta#eta<0.5)",
00469                                             8000,0.,8000.,
00470                                             "M_{jj} WideJets [GeV]",
00471                                             "Number of events"
00472                                             );
00473 
00474   m_MjjWide_deta_0p5_1p0 = bookH1withSumw2( "h1_MjjWide_deta_0p5_1p0",
00475                                             "M_{jj} WideJets (0.5<=#Delta#eta<1.0)",
00476                                             8000,0.,8000.,
00477                                             "M_{jj} WideJets [GeV]",
00478                                             "Number of events"
00479                                             );
00480 
00481   m_MjjWide_deta_1p0_1p5 = bookH1withSumw2( "h1_MjjWide_deta_1p0_1p5",
00482                                             "M_{jj} WideJets (1.0<=#Delta#eta<1.5)",
00483                                             8000,0.,8000.,
00484                                             "M_{jj} WideJets [GeV]",
00485                                             "Number of events"
00486                                             );
00487 
00488   m_MjjWide_deta_1p5_2p0 = bookH1withSumw2( "h1_MjjWide_deta_1p5_2p0",
00489                                             "M_{jj} WideJets (1.5<=#Delta#eta<2.0)",
00490                                             8000,0.,8000.,
00491                                             "M_{jj} WideJets [GeV]",
00492                                             "Number of events"
00493                                             );
00494 
00495   m_MjjWide_deta_2p0_2p5 = bookH1withSumw2( "h1_MjjWide_deta_2p0_2p5",
00496                                             "M_{jj} WideJets (2.0<=#Delta#eta<2.5)",
00497                                             8000,0.,8000.,
00498                                             "M_{jj} WideJets [GeV]",
00499                                             "Number of events"
00500                                             );
00501 
00502   m_MjjWide_deta_2p5_3p0 = bookH1withSumw2( "h1_MjjWide_deta_2p5_3p0",
00503                                             "M_{jj} WideJets (2.5<=#Delta#eta<3.0)",
00504                                             8000,0.,8000.,
00505                                             "M_{jj} WideJets [GeV]",
00506                                             "Number of events"
00507                                             );
00508 
00509   m_MjjWide_deta_3p0_inf = bookH1withSumw2( "h1_MjjWide_deta_3p0_inf",
00510                                             "M_{jj} WideJets (#Delta#eta>=3.0)",
00511                                             8000,0.,8000.,
00512                                             "M_{jj} WideJets [GeV]",
00513                                             "Number of events"
00514                                             );
00515 
00516   m_MjjWide_den_NOdeta = bookH1withSumw2( "h1_MjjWide_den_NOdeta",
00517                                           "HLT Efficiency Studies (no deta cut)",
00518                                           400,0.,2000.,
00519                                           "M_{jj} WideJets [GeV]",
00520                                           "Number of events"
00521                                           );
00522   
00523   m_MjjWide_num_NOdeta = bookH1withSumw2( "h1_MjjWide_num_NOdeta",
00524                                           "HLT Efficiency Studies (no deta cut)",
00525                                           400,0.,2000.,
00526                                           "M_{jj} WideJets [GeV]",
00527                                           "Number of events"
00528                                           );
00529 
00530   m_MjjWide_den_detaL4 = bookH1withSumw2( "h1_MjjWide_den_detaL4",
00531                                           "HLT Efficiency Studies (deta cut < 4.0)",
00532                                           400,0.,2000.,
00533                                           "M_{jj} WideJets [GeV]",
00534                                           "Number of events"
00535                                           );
00536   
00537   m_MjjWide_num_detaL4 = bookH1withSumw2( "h1_MjjWide_num_detaL4",
00538                                           "HLT Efficiency Studies (deta cut < 4.0)",
00539                                           400,0.,2000.,
00540                                           "M_{jj} WideJets [GeV]",
00541                                           "Number of events"
00542                                           );
00543 
00544   m_MjjWide_den_detaL3 = bookH1withSumw2( "h1_MjjWide_den_detaL3",
00545                                           "HLT Efficiency Studies (deta cut < 3.0)",
00546                                           400,0.,2000.,
00547                                           "M_{jj} WideJets [GeV]",
00548                                           "Number of events"
00549                                           );
00550   
00551   m_MjjWide_num_detaL3 = bookH1withSumw2( "h1_MjjWide_num_detaL3",
00552                                           "HLT Efficiency Studies (deta cut < 3.0)",
00553                                           400,0.,2000.,
00554                                           "M_{jj} WideJets [GeV]",
00555                                           "Number of events"
00556                                           );
00557 
00558   m_MjjWide_den_detaL2 = bookH1withSumw2( "h1_MjjWide_den_detaL2",
00559                                           "HLT Efficiency Studies (deta cut < 2.0)",
00560                                           400,0.,2000.,
00561                                           "M_{jj} WideJets [GeV]",
00562                                           "Number of events"
00563                                           );
00564   
00565   m_MjjWide_num_detaL2 = bookH1withSumw2( "h1_MjjWide_num_detaL2",
00566                                           "HLT Efficiency Studies (deta cut < 2.0)",
00567                                           400,0.,2000.,
00568                                           "M_{jj} WideJets [GeV]",
00569                                           "Number of events"
00570                                           );
00571 
00572   m_MjjWide_den = bookH1withSumw2( "h1_MjjWide_den",
00573                                    "HLT Efficiency Studies (default deta cut)",
00574                                    400,0.,2000.,
00575                                    "M_{jj} WideJets [GeV]",
00576                                    "Number of events"
00577                                    );
00578 
00579   m_MjjWide_num = bookH1withSumw2( "h1_MjjWide_num",
00580                                    "HLT Efficiency Studies (default deta cut)",
00581                                    400,0.,2000.,
00582                                    "M_{jj} WideJets [GeV]",
00583                                    "Number of events"
00584                                    );
00585 
00586   m_DetajjWide_finalSel = bookH1withSumw2( "h1_DetajjWide_finalSel",
00587                                            "#Delta#eta_{jj} WideJets (final selection)",
00588                                            100,0.,5.,
00589                                            "#Delta#eta_{jj} WideJets",
00590                                            "Number of events"
00591                                            );
00592 
00593   m_DetajjWide = bookH1withSumw2( "h1_DetajjWide",
00594                                   "#Delta#eta_{jj} WideJets (final selection except #Delta#eta cut)",
00595                                   100,0.,5.,
00596                                   "#Delta#eta_{jj} WideJets",
00597                                   "Number of events"
00598                                   );
00599 
00600   m_DphijjWide_finalSel = bookH1withSumw2( "h1_DphijjWide_finalSel",
00601                                            "#Delta#phi_{jj} WideJets (final selection)",
00602                                            100,0.,TMath::Pi()+0.0001,
00603                                            "#Delta#phi_{jj} WideJets [rad.]",
00604                                            "Number of events"
00605                                            );
00606 
00607 
00608   m_selJets_pt = bookH1withSumw2( "h1_selJets_pt",
00609                                   "Selected CaloJets",
00610                                   500,0.,5000.,
00611                                   "Jet Pt [GeV]",
00612                                   "Number of events"
00613                                   );
00614 
00615   m_selJets_eta = bookH1withSumw2( "h1_selJets_eta",
00616                                   "Selected CaloJets",
00617                                   100,-5.,5.,
00618                                   "#eta",
00619                                   "Number of events"
00620                                   );
00621 
00622   m_selJets_phi = bookH1withSumw2( "h1_selJets_phi",
00623                                   "Selected CaloJets",
00624                                    100,-TMath::Pi(),TMath::Pi(),
00625                                   "#phi (rad.)",
00626                                   "Number of events"
00627                                   );
00628 
00629   m_selJets_hadEnergyFraction = bookH1withSumw2( "h1_selJets_hadEnergyFraction",
00630                                                  "Selected CaloJets",
00631                                                  110,0.,1.1,
00632                                                  "HAD Energy Fraction",
00633                                                  "Number of events"
00634                                                  );
00635 
00636   m_selJets_emEnergyFraction = bookH1withSumw2( "h1_selJets_emEnergyFraction",
00637                                                  "Selected CaloJets",
00638                                                  110,0.,1.1,
00639                                                  "EM Energy Fraction",
00640                                                  "Number of events"
00641                                                  );
00642 
00643   m_selJets_towersArea = bookH1withSumw2( "h1_selJets_towersArea",
00644                                           "Selected CaloJets",
00645                                           200,0.,2.,
00646                                           "towers area",
00647                                           "Number of events"
00648                                           );
00649 
00650   m_metDiff = bookH1withSumw2( "h1_metDiff",
00651                                "Met - MetCleaned",
00652                                500,-1000.,1000.,
00653                                "met - metcleaned [GeV]",
00654                                "Number of events"
00655                                );
00656 
00657   m_metCases = bookH1withSumw2( "h1_metCases",
00658                                 "Met cases",
00659                                 3,0.,3.,
00660                                 "case",
00661                                 "Number of events"
00662                                 );
00663   m_metCases->getTH1()->GetXaxis()->SetBinLabel(1,"met , metclean");
00664   m_metCases->getTH1()->GetXaxis()->SetBinLabel(2,"met , !metclean");
00665   m_metCases->getTH1()->GetXaxis()->SetBinLabel(3,"!met , metclean");
00666 
00667 
00668   m_metCaseNoMetClean = bookH1withSumw2( "h1_metCaseNoMetClean",
00669                                          "Met - MetCleaned",
00670                                          1000,0.,2000.,
00671                                          "MET [GeV]",
00672                                          "Number of events"
00673                                          );  
00674 
00675   m_HT_inclusive = bookH1withSumw2( "h1_HT_inclusive",
00676                                     "HT (inclusive)",
00677                                     150,0.,15000.,
00678                                     "HT [GeV]",
00679                                     "Number of events"
00680                                     );  
00681 
00682   m_HT_finalSel = bookH1withSumw2( "h1_HT_finalSel",
00683                                    "HT (final selection)",
00684                                    150,0.,15000.,
00685                                    "HT [GeV]",
00686                                    "Number of events"
00687                                    );  
00688 
00689 
00690   // 2D histograms
00691   m_DetajjVsMjjWide = bookH2withSumw2("h2_DetajjVsMjjWide",
00692                                       "#Delta#eta_{jj} vs M_{jj} WideJets",
00693                                       8000,0.,8000.,
00694                                       100,0.,5.,
00695                                       "M_{jj} WideJets [GeV]",
00696                                       "#Delta#eta_{jj} WideJets");
00697 
00698   m_DetajjVsMjjWide_rebin = bookH2withSumw2("h2_DetajjVsMjjWide_rebin",
00699                                             "#Delta#eta_{jj} vs M_{jj} WideJets",
00700                                             400,0.,8000.,
00701                                             50,0.,5.,
00702                                             "M_{jj} WideJets [GeV]",
00703                                             "#Delta#eta_{jj} WideJets");
00704 
00705   m_metVSmetclean = bookH2withSumw2("h2_metVSmetclean",
00706                                     "MET clean vs MET",
00707                                     100,0.,2000.,
00708                                     100,0.,2000.,
00709                                     "MET [GeV]",
00710                                     "MET clean [GeV]");
00711 
00712 
00713 }
00714