CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/DQMOffline/Trigger/src/HLTInclusiveVBFSource.cc

Go to the documentation of this file.
00001 /*
00002 */
00003 #include "DQMOffline/Trigger/interface/HLTInclusiveVBFSource.h"
00004 
00005 #include "FWCore/Common/interface/TriggerNames.h"
00006 #include "FWCore/Framework/interface/EDAnalyzer.h"
00007 #include "FWCore/Framework/interface/Run.h"
00008 #include "FWCore/Framework/interface/MakerMacros.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012 #include "FWCore/ServiceRegistry/interface/Service.h"
00013 
00014 #include "DQMServices/Core/interface/MonitorElement.h"
00015 
00016 #include "DataFormats/Common/interface/Handle.h"
00017 #include "DataFormats/Common/interface/TriggerResults.h"
00018 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00019 #include "DataFormats/HLTReco/interface/TriggerObject.h"
00020 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00021 #include "DataFormats/Math/interface/deltaR.h"
00022 #include "DataFormats/Math/interface/deltaPhi.h"
00023 #include "DataFormats/VertexReco/interface/Vertex.h"
00024 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00025 
00026 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
00027 
00028 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
00029 
00030 #include "math.h"
00031 #include "TH1F.h"
00032 #include "TProfile.h"
00033 #include "TH2F.h"
00034 #include "TPRegexp.h"
00035 #include "TMath.h"
00036 
00037 using namespace edm;
00038 using namespace reco;
00039 using namespace std;
00040 
00041   
00042 HLTInclusiveVBFSource::HLTInclusiveVBFSource(const edm::ParameterSet& iConfig):
00043   isSetup_(false)
00044 {
00045   LogDebug("HLTInclusiveVBFSource") << "constructor....";
00046   nCount_ = 0;
00047   dbe = Service < DQMStore > ().operator->();
00048   if ( ! dbe ) {
00049     LogDebug("HLTInclusiveVBFSource") << "unabel to get DQMStore service?";
00050   }
00051   if (iConfig.getUntrackedParameter < bool > ("DQMStore", false)) {
00052     dbe->setVerbose(0);
00053   }
00054   
00055   dirname_             = iConfig.getUntrackedParameter("dirname",std::string("HLT/InclusiveVBF"));
00056   processname_         = iConfig.getParameter<std::string>("processname");
00057   triggerSummaryLabel_ = iConfig.getParameter<edm::InputTag>("triggerSummaryLabel");
00058   triggerResultsLabel_ = iConfig.getParameter<edm::InputTag>("triggerResultsLabel");
00059   //path_                = iConfig.getUntrackedParameter<std::vector<std::string> >("paths");
00060   //l1path_              = iConfig.getUntrackedParameter<std::vector<std::string> >("l1paths"); 
00061   
00062   debug_               = iConfig.getUntrackedParameter< bool >("debug", false);
00063   
00064   caloJetsTag_         = iConfig.getParameter<edm::InputTag>("CaloJetCollectionLabel");
00065   caloMETTag_          = iConfig.getParameter<edm::InputTag>("CaloMETCollectionLabel"); 
00066   pfJetsTag_           = iConfig.getParameter<edm::InputTag>("PFJetCollectionLabel");
00067   pfMetTag_            = iConfig.getParameter<edm::InputTag>("PFMETCollectionLabel");
00068   //jetID                = new reco::helper::JetIDHelper(iConfig.getParameter<ParameterSet>("JetIDParams"));
00069   
00070   minPtHigh_           = iConfig.getUntrackedParameter<double>("minPtHigh",40.);
00071   minPtLow_            = iConfig.getUntrackedParameter<double>("minPtLow",40.);
00072   minDeltaEta_         = iConfig.getUntrackedParameter<double>("minDeltaEta",3.5);
00073   deltaRMatch_         = iConfig.getUntrackedParameter<double>("deltaRMatch",0.1);
00074   minInvMass_          = iConfig.getUntrackedParameter<double>("minInvMass",1000.0);
00075   etaOpposite_         = iConfig.getUntrackedParameter<bool>("etaOpposite",true);
00076 
00077   check_mjj650_Pt35_DEta3p5 = false;
00078   check_mjj700_Pt35_DEta3p5 = false;
00079   check_mjj750_Pt35_DEta3p5 = false;
00080   check_mjj800_Pt35_DEta3p5 = false;
00081   check_mjj650_Pt40_DEta3p5 = false; 
00082   check_mjj700_Pt40_DEta3p5 = false;
00083   check_mjj750_Pt40_DEta3p5 = false; 
00084   check_mjj800_Pt40_DEta3p5 = false;
00085 }
00086 
00087 
00088 HLTInclusiveVBFSource::~HLTInclusiveVBFSource() {
00089   //
00090   // do anything here that needs to be done at desctruction time
00091   // (e.g. close files, deallocate resources etc.)
00092 }
00093 
00094 
00095 void
00096 HLTInclusiveVBFSource::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup){
00097   using namespace std;
00098   using namespace edm;
00099   using namespace trigger;
00100   using namespace reco;
00101 
00102   if(debug_) cout<<"DEBUG-0: Start to analyze"<<endl;
00103   
00104   //****************************************************
00105   // Get trigger information. 
00106   //****************************************************
00107   //
00108   //---------- triggerResults ----------
00109   iEvent.getByLabel(triggerResultsLabel_, triggerResults_);
00110   if(!triggerResults_.isValid()) {
00111     edm::InputTag triggerResultsLabelFU(triggerResultsLabel_.label(),triggerResultsLabel_.instance(), "FU");
00112     iEvent.getByLabel(triggerResultsLabelFU,triggerResults_);
00113     if(!triggerResults_.isValid()) {
00114       edm::LogInfo("FourVectorHLTOffline") << "TriggerResults not found, "
00115         "skipping event";
00116       return;
00117     }
00118   }
00119   //
00120   //int npath;
00121   if(&triggerResults_) {
00122     // Check how many HLT triggers are in triggerResults
00123     //npath = triggerResults_->size();
00124     triggerNames_ = iEvent.triggerNames(*triggerResults_);
00125   } 
00126   else {
00127     edm::LogInfo("CaloMETHLTOfflineSource") << "TriggerResults::HLT not found, "
00128       "automatically select events";
00129     return;
00130   }
00131   //
00132   //---------- triggerSummary ----------
00133   iEvent.getByLabel(triggerSummaryLabel_,triggerObj_);
00134   if(!triggerObj_.isValid()) {
00135     edm::InputTag triggerSummaryLabelFU(triggerSummaryLabel_.label(),triggerSummaryLabel_.instance(), "FU");
00136     iEvent.getByLabel(triggerSummaryLabelFU,triggerObj_);
00137     if(!triggerObj_.isValid()) {
00138       edm::LogInfo("FourVectorHLTOffline") << "TriggerEvent not found, "
00139         "skipping event";
00140       return;
00141     }
00142   }
00143   
00144   if(debug_) cout<<"DEBUG-1: Trigger information"<<endl;
00145   
00146   //****************************************************
00147   // Get AOD information
00148   //****************************************************
00149   //
00150   edm::Handle<edm::View<reco::PFMET> > metSrc;
00151   bool ValidPFMET_ = iEvent.getByLabel(pfMetTag_ , metSrc);
00152   if(!ValidPFMET_) return;
00153   
00154   edm::Handle<edm::View<reco::PFJet> > jetSrc;
00155   bool ValidPFJet_ = iEvent.getByLabel(pfJetsTag_ , jetSrc);
00156   if(!ValidPFJet_) return;
00157   
00158   if(!metSrc.isValid()) return;
00159   if(!jetSrc.isValid()) return;
00160   const edm::View<reco::PFMET> & mets = *metSrc;
00161   const edm::View<reco::PFJet> & jets = *jetSrc;
00162   if(jets.size()<=0)    return;
00163   if(mets.size()<=0)    return;
00164 
00165   if(debug_) cout<<"DEBUG-2: AOD Information"<<endl;
00166   
00167   //****************************************************
00168   // Variable setting
00169   //****************************************************
00170   //
00171   pathname   = "dummy";
00172   filtername = "dummy";
00173   
00174   //
00175   reco_ejet1                = 0.; 
00176   //reco_etjet1               = 0.;
00177   reco_pxjet1               = 0.;
00178   reco_pyjet1               = 0.;
00179   reco_pzjet1               = 0.;
00180   reco_ptjet1               = 0.;
00181   reco_etajet1              = 0.;
00182   reco_phijet1              = 0.;
00183   
00184   //
00185   reco_ejet2                = 0.;
00186   //reco_etjet2               = 0.;
00187   reco_pxjet2               = 0.;
00188   reco_pyjet2               = 0.;
00189   reco_pzjet2               = 0.;
00190   reco_ptjet2               = 0.;
00191   reco_etajet2              = 0.;
00192   reco_phijet2              = 0.;
00193   
00194   //
00195   hlt_ejet1                 = 0.;
00196   //hlt_etjet1                = 0.;
00197   hlt_pxjet1                = 0.;
00198   hlt_pyjet1                = 0.;
00199   hlt_pzjet1                = 0.;
00200   hlt_ptjet1                = 0.;
00201   hlt_etajet1               = 0.;
00202   hlt_phijet1               = 0.;
00203   
00204   //
00205   hlt_ejet2                 = 0.;
00206   //hlt_etjet2                = 0.;
00207   hlt_pxjet2                = 0.;
00208   hlt_pyjet2                = 0.;
00209   hlt_pzjet2                = 0.;
00210   hlt_ptjet2                = 0.;
00211   hlt_etajet2               = 0.;
00212   hlt_phijet2               = 0.;
00213   
00214   //
00215   checkOffline              = false;
00216   checkHLT                  = false;
00217   checkHLTIndex             = false;
00218   
00219   //
00220   dR_HLT_RECO_11            = 0.; 
00221   dR_HLT_RECO_22            = 0.;
00222   dR_HLT_RECO_12            = 0.; 
00223   dR_HLT_RECO_21            = 0.;
00224   
00225   //
00226   checkdR_sameOrder         = false;
00227   checkdR_crossOrder        = false;
00228   
00229   //
00230   reco_deltaetajet          = 0.;
00231   reco_deltaphijet          = 0.;
00232   reco_invmassjet           = 0.;
00233   hlt_deltaetajet           = 0.;
00234   hlt_deltaphijet           = 0.;
00235   hlt_invmassjet            = 0.;
00236   
00237   //****************************************************
00238   // Offline analysis
00239   //****************************************************
00240   //
00241   checkOffline  = false;
00242   for(unsigned int ijet1=0; ijet1<jets.size(); ijet1++){
00243     if(jets[ijet1].neutralHadronEnergyFraction()>0.99) continue;
00244     if(jets[ijet1].neutralEmEnergyFraction()>0.99) continue;
00245     for(unsigned int ijet2=ijet1+1; ijet2<jets.size(); ijet2++){ 
00246       if(jets[ijet2].neutralHadronEnergyFraction()>0.99) continue;
00247       if(jets[ijet2].neutralEmEnergyFraction()>0.99) continue;
00248       //
00249       reco_ejet1   = jets[ijet1].energy();
00250       //reco_etjet1  = jets[ijet1].et();
00251       reco_pxjet1  = jets[ijet1].momentum().X();
00252       reco_pyjet1  = jets[ijet1].momentum().Y();
00253       reco_pzjet1  = jets[ijet1].momentum().Z();
00254       reco_ptjet1  = jets[ijet1].pt();
00255       reco_etajet1 = jets[ijet1].eta();
00256       reco_phijet1 = jets[ijet1].phi(); 
00257       //
00258       reco_ejet2   = jets[ijet2].energy();
00259       //reco_etjet2  = jets[ijet2].et();
00260       reco_pxjet2  = jets[ijet2].momentum().X();
00261       reco_pyjet2  = jets[ijet2].momentum().Y();
00262       reco_pzjet2  = jets[ijet2].momentum().Z();
00263       reco_ptjet2  = jets[ijet2].pt();
00264       reco_etajet2 = jets[ijet2].eta();
00265       reco_phijet2 = jets[ijet2].phi();
00266       //
00267       reco_deltaetajet  = reco_etajet1 - reco_etajet2;
00268       reco_deltaphijet  = reco::deltaPhi(reco_phijet1,reco_phijet2);
00269       reco_invmassjet   = sqrt((reco_ejet1  + reco_ejet2)  * (reco_ejet1  + reco_ejet2)  - 
00270                                (reco_pxjet1 + reco_pxjet2) * (reco_pxjet1 + reco_pxjet2) - 
00271                                (reco_pyjet1 + reco_pyjet2) * (reco_pyjet1 + reco_pyjet2) - 
00272                                (reco_pzjet1 + reco_pzjet2) * (reco_pzjet1 + reco_pzjet2));
00273       
00274       //
00275       if(reco_ptjet1 < minPtHigh_) continue; 
00276       if(reco_ptjet2 < minPtLow_) continue;
00277       if(etaOpposite_ == true && reco_etajet1*reco_etajet2 > 0) continue;
00278       if(std::abs(reco_deltaetajet) < minDeltaEta_) continue;
00279       if(std::abs(reco_invmassjet)  < minInvMass_)  continue;
00280       
00281       //
00282       if(debug_) cout<<"DEBUG-3"<<endl;
00283       checkOffline  = true;
00284       break;
00285     }
00286     if(checkOffline == true) break;
00287   }
00288   if(checkOffline == false) return;
00289   
00290   //****************************************************
00291   // Trigger efficiency: Loop for all VBF paths
00292   //****************************************************
00293   //const unsigned int numberOfPaths(hltConfig_.size()); 
00294   const trigger::TriggerObjectCollection & toc(triggerObj_->getObjects()); 
00295   for(PathInfoCollection::iterator v = hltPathsAll_.begin(); v!= hltPathsAll_.end(); ++v ){
00296     checkHLT = false;
00297     checkHLTIndex = false;
00298     
00299     //
00300     v->getMEhisto_RECO_deltaEta_DiJet()->Fill(reco_deltaetajet);
00301     v->getMEhisto_RECO_deltaPhi_DiJet()->Fill(reco_deltaphijet);
00302     v->getMEhisto_RECO_invMass_DiJet()->Fill(reco_invmassjet);
00303     
00304     //
00305     if(debug_) cout<<"DEBUG-4-0: Path loops"<<endl;
00306     
00307     //
00308     if(isHLTPathAccepted(v->getPath())==false) continue;
00309     checkHLT = true;
00310     
00311     //
00312     if(debug_) cout<<"DEBUG-4-1: Path is accepted. Now we are looking for "<<v->getLabel()<<" module."<<endl;
00313     
00314     //
00315     edm::InputTag hltTag(v->getLabel(),"",processname_);
00316     const int hltIndex = triggerObj_->filterIndex(hltTag);
00317     if(hltIndex >= triggerObj_->sizeFilters()) continue;
00318     checkHLT = true;
00319     if(debug_) cout<<"DEBUG-4-2: HLT module "<<v->getLabel()<<" exists"<<endl;
00320     const trigger::Keys & khlt = triggerObj_->filterKeys(hltIndex);
00321     trigger::Keys::const_iterator kj = khlt.begin();
00322     for(; kj != khlt.end(); kj+=2){
00323       if(debug_) cout<<"DEBUG-5"<<endl;
00324       checkdR_sameOrder  = false;
00325       checkdR_crossOrder = false;           //
00326       hlt_ejet1   = toc[*kj].energy();
00327       //hlt_etjet1  = toc[*kj].et();
00328       hlt_pxjet1  = toc[*kj].px();
00329       hlt_pyjet1  = toc[*kj].py();
00330       hlt_pzjet1  = toc[*kj].pz();
00331       hlt_ptjet1  = toc[*kj].pt();
00332       hlt_etajet1 = toc[*kj].eta();
00333       hlt_phijet1 = toc[*kj].phi();  
00334       //
00335       hlt_ejet2   = toc[*(kj+1)].energy();
00336       //hlt_etjet2  = toc[*(kj+1)].et();
00337       hlt_pxjet2  = toc[*(kj+1)].px();
00338       hlt_pyjet2  = toc[*(kj+1)].py();
00339       hlt_pzjet2  = toc[*(kj+1)].pz();
00340       hlt_ptjet2  = toc[*(kj+1)].pt();
00341       hlt_etajet2 = toc[*(kj+1)].eta();
00342       hlt_phijet2 = toc[*(kj+1)].phi(); 
00343       //
00344       dR_HLT_RECO_11 = reco::deltaR(hlt_etajet1,hlt_phijet1,reco_etajet1,reco_phijet1);
00345       dR_HLT_RECO_22 = reco::deltaR(hlt_etajet2,hlt_phijet2,reco_etajet2,reco_phijet2);
00346       dR_HLT_RECO_12 = reco::deltaR(hlt_etajet1,hlt_phijet1,reco_etajet2,reco_phijet2);
00347       dR_HLT_RECO_21 = reco::deltaR(hlt_etajet2,hlt_phijet2,reco_etajet1,reco_phijet1);
00348       if(dR_HLT_RECO_11<deltaRMatch_ && dR_HLT_RECO_22<deltaRMatch_) checkdR_sameOrder = true;
00349       if(dR_HLT_RECO_12<deltaRMatch_ && dR_HLT_RECO_21<deltaRMatch_) checkdR_crossOrder = true;
00350       if(checkdR_sameOrder  == false && checkdR_crossOrder == false) continue; 
00351       checkHLTIndex = true;
00352       //
00353       if(debug_) cout<<"DEBUG-6: Match"<<endl;
00354       hlt_deltaetajet  = hlt_etajet1 - hlt_etajet2;
00355       hlt_deltaphijet  = reco::deltaPhi(hlt_phijet1,hlt_phijet2);
00356       hlt_invmassjet   = sqrt((hlt_ejet1  + hlt_ejet2)  * (hlt_ejet1  + hlt_ejet2)  - 
00357                               (hlt_pxjet1 + hlt_pxjet2) * (hlt_pxjet1 + hlt_pxjet2) - 
00358                               (hlt_pyjet1 + hlt_pyjet2) * (hlt_pyjet1 + hlt_pyjet2) - 
00359                               (hlt_pzjet1 + hlt_pzjet2) * (hlt_pzjet1 + hlt_pzjet2));
00360       v->getMEhisto_HLT_deltaEta_DiJet()->Fill(hlt_deltaetajet);
00361       v->getMEhisto_HLT_deltaPhi_DiJet()->Fill(hlt_deltaphijet);
00362       v->getMEhisto_HLT_invMass_DiJet()->Fill(hlt_invmassjet);
00363       //
00364       v->getMEhisto_RECO_deltaEta_DiJet_Match()->Fill(reco_deltaetajet);
00365       v->getMEhisto_RECO_deltaPhi_DiJet_Match()->Fill(reco_deltaphijet);
00366       v->getMEhisto_RECO_invMass_DiJet_Match()->Fill(reco_invmassjet);
00367       //
00368       v->getMEhisto_RECOHLT_deltaEta()->Fill(reco_deltaetajet,hlt_deltaetajet); 
00369       v->getMEhisto_RECOHLT_deltaPhi()->Fill(reco_deltaphijet,hlt_deltaphijet);
00370       v->getMEhisto_RECOHLT_invMass()->Fill(reco_invmassjet,hlt_invmassjet);
00371       //
00372       if(checkHLTIndex==true) break;
00373     }
00374     
00375     //****************************************************
00376     // Match information
00377     //****************************************************
00378     if(checkHLT==true && checkHLTIndex==true){
00379       if(debug_) cout<<"DEBUG-7: Match"<<endl;
00380       v->getMEhisto_NumberOfMatches()->Fill(1);
00381     }
00382     else{
00383       if(debug_) cout<<"DEBUG-8: Not match"<<endl;
00384       v->getMEhisto_NumberOfMatches()->Fill(0);
00385     }
00386   }
00387   
00388   
00389   //****************************************************
00390   // 
00391   //****************************************************
00392   for(PathInfoCollection::iterator v = hltPathsAll_.begin(); v!= hltPathsAll_.end(); ++v ){
00393     if(isHLTPathAccepted(v->getPath())==false) continue;
00394     if(debug_) cout<<"DEBUG-9: Loop for rate approximation: "<<v->getPath()<<endl;
00395     check_mjj650_Pt35_DEta3p5 = false;
00396     check_mjj700_Pt35_DEta3p5 = false;
00397     check_mjj750_Pt35_DEta3p5 = false;
00398     check_mjj800_Pt35_DEta3p5 = false;
00399     check_mjj650_Pt40_DEta3p5 = false; 
00400     check_mjj700_Pt40_DEta3p5 = false;
00401     check_mjj750_Pt40_DEta3p5 = false; 
00402     check_mjj800_Pt40_DEta3p5 = false;
00403     edm::InputTag hltTag(v->getLabel(),"",processname_);
00404     const int hltIndex = triggerObj_->filterIndex(hltTag);
00405     if(hltIndex >= triggerObj_->sizeFilters()) continue;
00406     const trigger::Keys & khlt = triggerObj_->filterKeys(hltIndex);
00407     trigger::Keys::const_iterator kj = khlt.begin();
00408     for(; kj != khlt.end(); kj+=2){
00409       checkdR_sameOrder  = false;
00410       checkdR_crossOrder = false;
00411       //
00412       hlt_ejet1   = toc[*kj].energy();
00413       //hlt_etjet1  = toc[*kj].et();
00414       hlt_pxjet1  = toc[*kj].px();
00415       hlt_pyjet1  = toc[*kj].py();
00416       hlt_pzjet1  = toc[*kj].pz();
00417       hlt_ptjet1  = toc[*kj].pt();
00418       hlt_etajet1 = toc[*kj].eta();
00419       hlt_phijet1 = toc[*kj].phi();  
00420       //
00421       hlt_ejet2   = toc[*(kj+1)].energy();
00422       //hlt_etjet2  = toc[*(kj+1)].et();
00423       hlt_pxjet2  = toc[*(kj+1)].px();
00424       hlt_pyjet2  = toc[*(kj+1)].py();
00425       hlt_pzjet2  = toc[*(kj+1)].pz();
00426       hlt_ptjet2  = toc[*(kj+1)].pt();
00427       hlt_etajet2 = toc[*(kj+1)].eta();
00428       hlt_phijet2 = toc[*(kj+1)].phi(); 
00429       //
00430       hlt_deltaetajet  = hlt_etajet1 - hlt_etajet2;
00431       hlt_deltaphijet  = reco::deltaPhi(hlt_phijet1,hlt_phijet2);
00432       hlt_invmassjet   = sqrt((hlt_ejet1  + hlt_ejet2)  * (hlt_ejet1  + hlt_ejet2)  - 
00433                               (hlt_pxjet1 + hlt_pxjet2) * (hlt_pxjet1 + hlt_pxjet2) - 
00434                               (hlt_pyjet1 + hlt_pyjet2) * (hlt_pyjet1 + hlt_pyjet2) - 
00435                               (hlt_pzjet1 + hlt_pzjet2) * (hlt_pzjet1 + hlt_pzjet2));
00436       //
00437       if(check_mjj650_Pt35_DEta3p5==false && hlt_ptjet1>35. && hlt_ptjet2>=35. && hlt_invmassjet>650 && std::abs(hlt_deltaetajet)>3.5){
00438         check_mjj650_Pt35_DEta3p5=true;
00439       }
00440       if(check_mjj700_Pt35_DEta3p5==false && hlt_ptjet1>35. && hlt_ptjet2>=35. && hlt_invmassjet>700 && std::abs(hlt_deltaetajet)>3.5){
00441         check_mjj700_Pt35_DEta3p5=true;
00442       } 
00443       if(check_mjj750_Pt35_DEta3p5==false && hlt_ptjet1>35. && hlt_ptjet2>=35. && hlt_invmassjet>750 && std::abs(hlt_deltaetajet)>3.5){
00444         check_mjj750_Pt35_DEta3p5=true;
00445       }
00446       if(check_mjj800_Pt35_DEta3p5==false && hlt_ptjet1>35. && hlt_ptjet2>=35. && hlt_invmassjet>800 && std::abs(hlt_deltaetajet)>3.5){
00447         check_mjj800_Pt35_DEta3p5=true;
00448       }
00449       if(check_mjj650_Pt40_DEta3p5==false && hlt_ptjet1>40. && hlt_ptjet2>=40. && hlt_invmassjet>650 && std::abs(hlt_deltaetajet)>3.5){
00450         check_mjj650_Pt40_DEta3p5=true;
00451       }
00452       if(check_mjj700_Pt40_DEta3p5==false && hlt_ptjet1>40. && hlt_ptjet2>=40. && hlt_invmassjet>700 && std::abs(hlt_deltaetajet)>3.5){
00453         check_mjj700_Pt40_DEta3p5=true;
00454       } 
00455       if(check_mjj750_Pt40_DEta3p5==false && hlt_ptjet1>40. && hlt_ptjet2>=40. && hlt_invmassjet>750 && std::abs(hlt_deltaetajet)>3.5){
00456         check_mjj750_Pt40_DEta3p5=true;
00457       }
00458       if(check_mjj800_Pt40_DEta3p5==false && hlt_ptjet1>40. && hlt_ptjet2>=40. && hlt_invmassjet>800 && std::abs(hlt_deltaetajet)>3.5){
00459         check_mjj800_Pt40_DEta3p5=true;
00460       }
00461     }
00462     if(check_mjj650_Pt35_DEta3p5==true) v->getMEhisto_NumberOfEvents()->Fill(0);
00463     if(check_mjj700_Pt35_DEta3p5==true) v->getMEhisto_NumberOfEvents()->Fill(1);
00464     if(check_mjj750_Pt35_DEta3p5==true) v->getMEhisto_NumberOfEvents()->Fill(2);
00465     if(check_mjj800_Pt35_DEta3p5==true) v->getMEhisto_NumberOfEvents()->Fill(3);
00466     if(check_mjj650_Pt40_DEta3p5==true) v->getMEhisto_NumberOfEvents()->Fill(4);
00467     if(check_mjj700_Pt40_DEta3p5==true) v->getMEhisto_NumberOfEvents()->Fill(5);
00468     if(check_mjj750_Pt40_DEta3p5==true) v->getMEhisto_NumberOfEvents()->Fill(6);
00469     if(check_mjj800_Pt40_DEta3p5==true) v->getMEhisto_NumberOfEvents()->Fill(7);
00470   }
00471 }
00472 
00473 
00474 // -- method called once each job just before starting event loop  --------
00475 void 
00476 HLTInclusiveVBFSource::beginJob(){
00477 }
00478 
00479 // BeginRun
00480 void 
00481 HLTInclusiveVBFSource::beginRun(const edm::Run& run, const edm::EventSetup& c){
00482   if(!isSetup_){
00483     DQMStore *dbe = 0;
00484     dbe = Service<DQMStore>().operator->();
00485     if (dbe) {
00486       dbe->setCurrentFolder(dirname_);
00487       dbe->rmdir(dirname_);
00488     }
00489     if (dbe) {
00490       dbe->setCurrentFolder(dirname_);
00491     }
00492     
00493     //--- htlConfig_
00494     bool changed(true);
00495     if (!hltConfig_.init(run, c, processname_, changed)) {
00496       LogDebug("HLTInclusiveVBFSource") << "HLTConfigProvider failed to initialize.";
00497     }
00498     
00499     const unsigned int numberOfPaths(hltConfig_.size());
00500     for(unsigned int i=0; i!=numberOfPaths; ++i){
00501       bool numFound = false;
00502       pathname      = hltConfig_.triggerName(i);
00503       filtername    = "dummy";
00504       unsigned int usedPrescale = 1;
00505       unsigned int objectType = 0;
00506       std::string triggerType = "";
00507       
00508       if(pathname.find("HLT_Di") == std::string::npos) continue;
00509       if(pathname.find("Jet") == std::string::npos) continue;
00510       if(pathname.find("MJJ") == std::string::npos) continue;
00511       if(pathname.find("VBF_v") == std::string::npos) continue;
00512       
00513       if(debug_){
00514         cout<<" - Startup:Path = "<<pathname<<endl;
00515       //cout<<" - Startup:PS = "<<hltConfig_.prescaleSize()<<endl;
00516       }
00517       
00518       triggerType = "DiJet_Trigger";
00519       objectType = trigger::TriggerJet;
00520       
00521       // Checking if the trigger exist in HLT table or not
00522       for (unsigned int i=0; i!=numberOfPaths; ++i) {
00523         std::string HLTname = hltConfig_.triggerName(i);
00524         if(HLTname == pathname)numFound = true;
00525       }
00526       
00527       if(numFound==false) continue;
00528       std::vector<std::string> numpathmodules = hltConfig_.moduleLabels(pathname);
00529       std::vector<std::string>::iterator numpathmodule = numpathmodules.begin();
00530       for(; numpathmodule!= numpathmodules.end(); ++numpathmodule){
00531         edm::InputTag testTag(*numpathmodule,"",processname_);
00532         if (hltConfig_.moduleType(*numpathmodule) == "HLTCaloJetVBFFilter"
00533             || hltConfig_.moduleType(*numpathmodule) == "HLTPFJetVBFFilter")
00534           {
00535             filtername = *numpathmodule;
00536             if(debug_) cout<<" - Startup:Module = "<<hltConfig_.moduleType(*numpathmodule)<<", FilterName = "<<filtername<<endl;
00537           }
00538         
00539       }
00540       if(debug_) cout<<" - Startup:Final filter = "<<filtername<<endl;
00541       
00542       if(objectType == 0 || numFound==false) continue;
00543       //if(debug_){
00544       //cout<<"Pathname = "<<pathname
00545       //    <<", Filtername = "<<filtername
00546       //    <<", ObjectType = "<<objectType<<endl;
00547       //}
00548       hltPathsAll_.push_back(PathInfo(usedPrescale, pathname, filtername, processname_, objectType, triggerType)); 
00549     }//Loop over paths
00550     
00551     //if(debug_) cout<<"== end hltPathsEff_.push_back ======" << endl;
00552     
00553     std::string dirName = dirname_ + "/MonitorInclusiveVBFTrigger/";
00554     for(PathInfoCollection::iterator v = hltPathsAll_.begin(); v!= hltPathsAll_.end(); ++v ){
00555       if(debug_) cout<<"Storing: "<<v->getPath()<<", Prescale = "<<v->getprescaleUsed()<<endl;
00556       //if(v->getprescaleUsed()!=1) continue;
00557       
00558       std::string subdirName = dirName + v->getPath();
00559       std::string trigPath = "("+v->getPath()+")";
00560       dbe->setCurrentFolder(subdirName);  
00561       
00562       MonitorElement*  RECO_deltaEta_DiJet;
00563       MonitorElement*  RECO_deltaPhi_DiJet;
00564       MonitorElement*  RECO_invMass_DiJet;
00565       MonitorElement*  HLT_deltaEta_DiJet;
00566       MonitorElement*  HLT_deltaPhi_DiJet;
00567       MonitorElement*  HLT_invMass_DiJet;
00568       MonitorElement*  RECO_deltaEta_DiJet_Match;
00569       MonitorElement*  RECO_deltaPhi_DiJet_Match;
00570       MonitorElement*  RECO_invMass_DiJet_Match;
00571       MonitorElement*  RECOHLT_deltaEta;
00572       MonitorElement*  RECOHLT_deltaPhi;
00573       MonitorElement*  RECOHLT_invMass;
00574       MonitorElement*  NumberOfMatches;
00575       MonitorElement*  NumberOfEvents;
00576       
00577       //dummy                     = dbe->bookFloat("dummy");
00578       RECO_deltaEta_DiJet       = dbe->bookFloat("RECO_deltaEta_DiJet");
00579       RECO_deltaPhi_DiJet       = dbe->bookFloat("RECO_deltaPhi_DiJet");
00580       RECO_invMass_DiJet        = dbe->bookFloat("RECO_invMass_DiJet");
00581       HLT_deltaEta_DiJet        = dbe->bookFloat("HLT_deltaEta_DiJet");
00582       HLT_deltaPhi_DiJet        = dbe->bookFloat("HLT_deltaPhi_DiJet ");
00583       HLT_invMass_DiJet         = dbe->bookFloat("HLT_invMass_DiJet");
00584       RECO_deltaEta_DiJet_Match = dbe->bookFloat("RECO_deltaEta_DiJet_Match");
00585       RECO_deltaPhi_DiJet_Match = dbe->bookFloat("RECO_deltaPhi_DiJet_Match");
00586       RECO_invMass_DiJet_Match  = dbe->bookFloat("RECO_invMass_DiJet_Match");
00587       RECOHLT_deltaEta          = dbe->bookFloat("RECOHLT_deltaEta");
00588       RECOHLT_deltaPhi          = dbe->bookFloat("RECOHLT_deltaPhi ");
00589       RECOHLT_invMass           = dbe->bookFloat("RECOHLT_invMass");
00590       NumberOfMatches           = dbe->bookFloat("NumberOfMatches");
00591       NumberOfEvents            = dbe->bookFloat("NumberOfEvents");
00592       
00593       std::string labelname("ME");
00594       std::string histoname(labelname+"");
00595       std::string title(labelname+"");
00596         
00597       //RECO_deltaEta_DiJet
00598       histoname     = labelname+"_RECO_deltaEta_DiJet";
00599       title         = labelname+"_RECO_deltaEta_DiJet "+trigPath;
00600       RECO_deltaEta_DiJet = dbe->book1D(histoname.c_str(),title.c_str(),50,-10.,10.);
00601       RECO_deltaEta_DiJet->getTH1F();
00602 
00603       //RECO_deltaPhi_DiJet
00604       histoname     = labelname+"_RECO_deltaPhi_DiJet";
00605       title         = labelname+"_RECO_deltaPhi_DiJet "+trigPath;
00606       RECO_deltaPhi_DiJet = dbe->book1D(histoname.c_str(),title.c_str(),35,-3.5,3.5);
00607       RECO_deltaPhi_DiJet->getTH1F();
00608 
00609       //RECO_invMass_DiJet
00610       histoname     = labelname+"_RECO_invMass_DiJet";
00611       title         = labelname+"_RECO_invMass_DiJet "+trigPath;
00612       RECO_invMass_DiJet = dbe->book1D(histoname.c_str(),title.c_str(),100,500.,2000.);
00613       RECO_invMass_DiJet->getTH1F(); 
00614 
00615       //HLT_deltaEta_DiJet
00616       histoname     = labelname+"_HLT_deltaEta_DiJet";
00617       title         = labelname+"_HLT_deltaEta_DiJet "+trigPath;
00618       HLT_deltaEta_DiJet = dbe->book1D(histoname.c_str(),title.c_str(),50,-10.,10.);
00619       HLT_deltaEta_DiJet->getTH1F();
00620 
00621       //HLT_deltaPhi_DiJet
00622       histoname     = labelname+"_HLT_deltaPhi_DiJet";
00623       title         = labelname+"_HLT_deltaPhi_DiJet "+trigPath;
00624       HLT_deltaPhi_DiJet = dbe->book1D(histoname.c_str(),title.c_str(),35,-3.5,3.5);
00625       HLT_deltaPhi_DiJet->getTH1F();
00626 
00627       //HLT_invMass_DiJet
00628       histoname     = labelname+"_HLT_invMass_DiJet";
00629       title         = labelname+"_HLT_invMass_DiJet "+trigPath;
00630       HLT_invMass_DiJet = dbe->book1D(histoname.c_str(),title.c_str(),100,500.,2000.);
00631       HLT_invMass_DiJet->getTH1F();
00632 
00633       //RECO_deltaEta_DiJet_Match
00634       histoname     = labelname+"_RECO_deltaEta_DiJet_Match";
00635       title         = labelname+"_RECO_deltaEta_DiJet_Match "+trigPath;
00636       RECO_deltaEta_DiJet_Match = dbe->book1D(histoname.c_str(),title.c_str(),50,-10.,10.);
00637       RECO_deltaEta_DiJet_Match->getTH1F();
00638 
00639       //RECO_deltaPhi_DiJet_Match
00640       histoname     = labelname+"_RECO_deltaPhi_DiJet_Match";
00641       title         = labelname+"_RECO_deltaPhi_DiJet_Match "+trigPath;
00642       RECO_deltaPhi_DiJet_Match = dbe->book1D(histoname.c_str(),title.c_str(),35,-3.5,3.5);
00643       RECO_deltaPhi_DiJet_Match->getTH1F();
00644 
00645       //RECO_invMass_DiJet_Match
00646       histoname     = labelname+"_RECO_invMass_DiJet_Match";
00647       title         = labelname+"_RECO_invMass_DiJet_Match "+trigPath;
00648       RECO_invMass_DiJet_Match = dbe->book1D(histoname.c_str(),title.c_str(),100,500.,2000.);
00649       RECO_invMass_DiJet_Match->getTH1F(); 
00650 
00651       //RECOHLT_deltaEta
00652       histoname     = labelname+"_RECOHLT_deltaEta";
00653       title         = labelname+"_RECOHLT_deltaEta "+trigPath;
00654       RECOHLT_deltaEta = dbe->book2D(histoname.c_str(),title.c_str(),50,-10.,10.,50,-10.,10.);
00655       RECOHLT_deltaEta->getTH2F();
00656 
00657       //RECOHLT_deltaPhi
00658       histoname     = labelname+"_RECOHLT_deltaPhi";
00659       title         = labelname+"_RECOHLT_deltaPhi "+trigPath;
00660       RECOHLT_deltaPhi = dbe->book2D(histoname.c_str(),title.c_str(),35,-3.5,3.5,35,-3.5,3.5);
00661       RECOHLT_deltaPhi->getTH2F();
00662 
00663       //RECOHLT_invMass
00664       histoname     = labelname+"_RECOHLT_invMass";
00665       title         = labelname+"_RECOHLT_invMass "+trigPath;
00666       RECOHLT_invMass = dbe->book2D(histoname.c_str(),title.c_str(),100,500.,2000.,100,500.,2000.);
00667       RECOHLT_invMass->getTH2F(); 
00668       
00669       //NumberOfMatches 
00670       histoname     = labelname+"_NumberOfMatches ";
00671       title         = labelname+"_NumberOfMatches  "+trigPath;
00672       NumberOfMatches = dbe->book1D(histoname.c_str(),title.c_str(),2,0.,2.);
00673       NumberOfMatches->getTH1F();
00674 
00675       //NumberOfEvents
00676       histoname     = labelname+"_NumberOfEvents";
00677       title         = labelname+"_NumberOfEvents "+trigPath;
00678       NumberOfEvents = dbe->book1D(histoname.c_str(),title.c_str(),10,0.,10.);
00679       NumberOfEvents->getTH1F();
00680       
00681       //} 
00682       v->setHistos(
00683                    RECO_deltaEta_DiJet,
00684                    RECO_deltaPhi_DiJet,
00685                    RECO_invMass_DiJet,
00686                    HLT_deltaEta_DiJet,
00687                    HLT_deltaPhi_DiJet,
00688                    HLT_invMass_DiJet,
00689                    RECO_deltaEta_DiJet_Match,
00690                    RECO_deltaPhi_DiJet_Match,
00691                    RECO_invMass_DiJet_Match,
00692                    RECOHLT_deltaEta,
00693                    RECOHLT_deltaPhi,
00694                    RECOHLT_invMass,
00695                    NumberOfMatches,
00696                    NumberOfEvents
00697                    );
00698       //break;//We need only the first unprescale paths
00699     }
00700   }
00701 }
00702 
00703 //--------------------------------------------------------
00704 void HLTInclusiveVBFSource::beginLuminosityBlock(const LuminosityBlock& lumiSeg, 
00705                                                   const EventSetup& context) {
00706 }
00707 
00708 //--------------------------------------------------------
00709 void 
00710 HLTInclusiveVBFSource::endLuminosityBlock(const LuminosityBlock& lumiSeg, 
00711                                                 const EventSetup& context) {
00712 }
00713 
00714 // - method called once each job just after ending the event loop  ------------
00715 void 
00716 HLTInclusiveVBFSource::endJob() {
00717   //delete jetID;
00718 }
00719 
00721 void HLTInclusiveVBFSource::endRun(const edm::Run& run, const edm::EventSetup& c){
00722   //if (debug_) std::cout << "endRun, run " << run.id() << std::endl;
00723 }
00724 
00725 bool HLTInclusiveVBFSource::isBarrel(double eta){
00726   bool output = false;
00727   if (fabs(eta)<=1.3) output=true;
00728   return output;
00729 }
00730 
00731 bool HLTInclusiveVBFSource::isEndCap(double eta){
00732   bool output = false;
00733   if (fabs(eta)<=3.0 && fabs(eta)>1.3) output=true;
00734   return output;
00735 }
00736 
00737 bool HLTInclusiveVBFSource::isForward(double eta){
00738   bool output = false;
00739   if (fabs(eta)>3.0) output=true;
00740   return output;
00741 }
00742 
00743 bool HLTInclusiveVBFSource::validPathHLT(std::string pathname){
00744   // hltConfig_ has to be defined first before calling this method
00745   bool output=false;
00746   for (unsigned int j=0; j!=hltConfig_.size(); ++j) {
00747     if (hltConfig_.triggerName(j) == pathname )
00748       output=true;
00749   }
00750   return output;
00751 }
00752 
00753 bool HLTInclusiveVBFSource::isHLTPathAccepted(std::string pathName){
00754   // triggerResults_, triggerNames_ has to be defined first before calling this method
00755   bool output=false;
00756   if(&triggerResults_) {
00757     unsigned index = triggerNames_.triggerIndex(pathName);
00758     //std::cout<<" -index = "<<index<<endl;
00759     if(index < triggerNames_.size() && triggerResults_->accept(index)) output = true;
00760   }
00761   return output;
00762 }
00763 
00764 bool HLTInclusiveVBFSource::isTriggerObjectFound(std::string objectName){
00765   // processname_, triggerObj_ has to be defined before calling this method
00766   bool output=false;
00767   edm::InputTag testTag(objectName,"",processname_);
00768   const int index = triggerObj_->filterIndex(testTag);
00769   if ( index >= triggerObj_->sizeFilters() ) {    
00770     edm::LogInfo("HLTInclusiveVBFSource") << "no index "<< index << " of that name ";
00771   } else {       
00772     const trigger::Keys & k = triggerObj_->filterKeys(index);
00773     if (k.size()) output=true;
00774   }
00775   return output;
00776 }
00777 
00778