CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/DQMOffline/Trigger/src/TopHLTDiMuonDQM.cc

Go to the documentation of this file.
00001 /*
00002  *  $Date: 2010/08/13 09:11:38 $
00003  *  $Revision: 1.8 $
00004  *  \author M. Marienfeld - DESY Hamburg
00005  */
00006 
00007 #include "DQMOffline/Trigger/interface/TopHLTDiMuonDQM.h"
00008 #include "FWCore/Common/interface/TriggerNames.h"
00009 
00010 using namespace std;
00011 using namespace edm;
00012 using namespace trigger;
00013 
00014 
00015 TopHLTDiMuonDQM::TopHLTDiMuonDQM( const edm::ParameterSet& ps ) {
00016 
00017   monitorName_ = ps.getParameter<string>("monitorName");
00018 
00019   triggerResults_ = ps.getParameter<edm::InputTag>("TriggerResults");
00020   triggerEvent_   = ps.getParameter<edm::InputTag>("TriggerEvent");
00021   triggerFilter_  = ps.getParameter<edm::InputTag>("TriggerFilter");
00022 
00023   hltPaths_L1_   = ps.getParameter<vector<string> >("hltPaths_L1");
00024   hltPaths_L3_   = ps.getParameter<vector<string> >("hltPaths_L3");
00025   hltPaths_sig_  = ps.getParameter<vector<string> >("hltPaths_sig");
00026   hltPaths_trig_ = ps.getParameter<vector<string> >("hltPaths_trig");
00027 
00028   vertex_       = ps.getParameter<edm::InputTag>("vertexCollection");
00029   vertex_X_cut_ = ps.getParameter<double>("vertex_X_cut");
00030   vertex_Y_cut_ = ps.getParameter<double>("vertex_Y_cut");
00031   vertex_Z_cut_ = ps.getParameter<double>("vertex_Z_cut");
00032 
00033   muons_        = ps.getParameter<edm::InputTag>("muonCollection");
00034   muon_pT_cut_  = ps.getParameter<double>("muon_pT_cut");
00035   muon_eta_cut_ = ps.getParameter<double>("muon_eta_cut");
00036   muon_iso_cut_ = ps.getParameter<double>("muon_iso_cut");
00037 
00038   MassWindow_up_   = ps.getParameter<double>("MassWindow_up");
00039   MassWindow_down_ = ps.getParameter<double>("MassWindow_down");
00040 
00041 }
00042 
00043 
00044 TopHLTDiMuonDQM::~TopHLTDiMuonDQM() {
00045 
00046 }
00047 
00048 
00049 void TopHLTDiMuonDQM::beginJob() {
00050 
00051   dbe_ = Service<DQMStore>().operator->();
00052 
00053   if( dbe_ ) {
00054 
00055     dbe_->setCurrentFolder(monitorName_);
00056 
00057     Trigs = dbe_->book1D("Trigs", "Fired triggers", 15, 0., 15.);
00058 
00059     TriggerEfficiencies = dbe_->book1D("TriggerEfficiencies", "HL Trigger Efficiencies", 10, 0., 10.);
00060 
00061     TriggerEfficiencies->setTitle("HL Trigger Efficiencies #epsilon_{signal} = #frac{[signal] && [control]}{[control]}");
00062 
00063     TriggerEfficiencies_sig  = dbe_->book1D("TriggerEfficiencies_sig",  "HL Trigger Signal && Control Counts", 10, 0., 10.);
00064     TriggerEfficiencies_trig = dbe_->book1D("TriggerEfficiencies_trig", "HL Trigger Control Counts",           10, 0., 10.);
00065 
00066     const int nbins_Pt = 5;
00067 
00068     float bins_Pt[nbins_Pt+1] = { 0., 5., 10., 20., 50., 100. };
00069 
00070     MuonEfficiency_pT      = dbe_->book1D("MuonEfficiency_pT",      "Muon Efficiency P_{T}",           nbins_Pt, &bins_Pt[0]);
00071     MuonEfficiency_pT_sig  = dbe_->book1D("MuonEfficiency_pT_sig",  "P^{#mu}_{T} (signal triggered)",  nbins_Pt, &bins_Pt[0]);
00072     MuonEfficiency_pT_trig = dbe_->book1D("MuonEfficiency_pT_trig", "P^{#mu}_{T} (control triggered)", nbins_Pt, &bins_Pt[0]);
00073 
00074     const int nbins_eta = 7;
00075 
00076     float bins_eta[nbins_eta+1] = { -3.5, -2.5, -1.5, -0.5, 0.5, 1.5, 2.5, 3.5 };
00077 
00078     MuonEfficiency_eta      = dbe_->book1D("MuonEfficiency_eta",      "Muon Efficiency  #eta",           nbins_eta, &bins_eta[0]);
00079     MuonEfficiency_eta_sig  = dbe_->book1D("MuonEfficiency_eta_sig",  "#eta_{muon} (signal triggered)",  nbins_eta, &bins_eta[0]);
00080     MuonEfficiency_eta_trig = dbe_->book1D("MuonEfficiency_eta_trig", "#eta_{muon} (control triggered)", nbins_eta, &bins_eta[0]);
00081 
00082     const int nbins_phi = 9;
00083 
00084     float bins_phi[nbins_phi+1] = { -3.5, -3.2, -2.6, -1.56, -0.52, 0.52, 1.56, 2.6, 3.2, 3.5 };
00085 
00086     MuonEfficiency_phi      = dbe_->book1D("MuonEfficiency_phi",      "Muon Efficiency  #phi",           nbins_phi, &bins_phi[0]);
00087     MuonEfficiency_phi_sig  = dbe_->book1D("MuonEfficiency_phi_sig",  "#phi_{muon} (signal triggered)",  nbins_phi, &bins_phi[0]);
00088     MuonEfficiency_phi_trig = dbe_->book1D("MuonEfficiency_phi_trig", "#phi_{muon} (control triggered)", nbins_phi, &bins_phi[0]);
00089 
00090     const int N_TriggerPaths = hltPaths_L3_.size();
00091     const int N_SignalPaths  = hltPaths_sig_.size();
00092 
00093     for( int i = 0; i < N_TriggerPaths; i++ ) {
00094 
00095       const string &label = hltPaths_L3_[i];
00096       Trigs->setBinLabel( i+1,label.c_str() );
00097 
00098     }
00099 
00100     for( int j = 0; j < N_SignalPaths; ++j ) {
00101 
00102       const string &label_eff = "#frac{["+hltPaths_sig_[j]+"]}{vs. ["+hltPaths_trig_[j]+"]}";
00103       const string &label_sig = hltPaths_sig_[j]+"\n && "+hltPaths_trig_[j];
00104       TriggerEfficiencies->setBinLabel(     j+1, label_eff.c_str() );
00105       TriggerEfficiencies_sig->setBinLabel( j+1, label_sig.c_str() );
00106 
00107     }
00108 
00109     NMuons        = dbe_->book1D("Nmuons",        "Number of muons",             20,   0.,  10.);
00110     NMuons_iso    = dbe_->book1D("Nmuons_Iso",    "Number of isolated muons",    20,   0.,  10.);
00111     NMuons_charge = dbe_->book1D("Nmuons_Charge", "N_{muons} * Q(#mu)",          19, -10.,  10.);
00112     NTracks       = dbe_->book1D("Ntracks",       "Number of tracks",            50,   0.,  50.);
00113     VxVy_muons    = dbe_->book2D("VxVy_muons",    "Vertex x-y-positon (global)", 40,  -1.,   1., 40 , -1., 1.);
00114     Vz_muons      = dbe_->book1D("Vz_muons",      "Vertex z-positon (global)",   40, -20.,  20.);
00115     PtMuons       = dbe_->book1D("PtMuon",        "P^{#mu}_{T}",                 50,   0., 100.);
00116     EtaMuons      = dbe_->book1D("EtaMuon",       "#eta_{muon}",                 50,  -5.,   5.);
00117     PhiMuons      = dbe_->book1D("PhiMuon",       "#phi_{muon}",                 40,  -4.,   4.);
00118 
00119     DeltaEtaMuonsRC   = dbe_->book1D("DeltaEtaMuonsRC",    "#Delta #eta of muon pair (RC)", 50,  -5.,   5.);
00120     DeltaPhiMuonsRC   = dbe_->book1D("DeltaPhiMuonsRC",    "#Delta #phi of muon pair (RC)", 50,  -5.,   5.);
00121     DeltaEtaMuonsWC   = dbe_->book1D("DeltaEtaMuonsWC",    "#Delta #eta of muon pair (WC)", 50,  -5.,   5.);
00122     DeltaPhiMuonsWC   = dbe_->book1D("DeltaPhiMuonsWC",    "#Delta #phi of muon pair (WC)", 50,  -5.,   5.);
00123     CombRelIso03      = dbe_->book1D("MuIso_CombRelIso03", "Muon CombRelIso dR=03",         50,   0.,   1.);
00124     PixelHits_muons   = dbe_->book1D("NPixelHits_muons",   "Number of pixel hits",          50,   0.,  50.);
00125     TrackerHits_muons = dbe_->book1D("NTrackerHits_muons", "Number of hits in the tracker", 50,   0.,  50.);
00126 
00127     DeltaR_Trig   = dbe_->book1D("DeltaRTrigger", "#Delta R of trigger muon pair",       50, 0., 5.);
00128     DeltaR_Reco   = dbe_->book1D("DeltaRReco",    "#Delta R of RECO muon pair",          50, 0., 5.);
00129     DeltaR_Match  = dbe_->book1D("DeltaRMatch",   "#Delta R of matched muon pair",       50, 0., 5.);
00130     Trigger_Match = dbe_->book1D("TriggerMatch",  "Number of Trigger-RECO assignements",  6, 0., 6.);
00131 
00132     const int nbins_Pt_Log = 15;
00133 
00134     double logmin = 0.;
00135     double logmax = 3.;  // 10^(3.)=1000
00136 
00137     float bins_Pt_Log[nbins_Pt_Log+1];
00138 
00139     for (int i = 0; i <= nbins_Pt_Log; i++) {
00140       double log = logmin + (logmax-logmin)*i/nbins_Pt_Log;
00141       bins_Pt_Log[i] = std::pow(10.0, log);
00142     }
00143 
00144     PtMuons_LOGX  = dbe_->book1D("Pt_muon_LOGX", "P^{#mu}_{T}", nbins_Pt_Log, &bins_Pt_Log[0]);
00145 
00146     MuonEfficiency_pT_LOGX      = dbe_->book1D("MuonEfficiency_pT_LOGX",      "Muon Efficiency P_{T}",           nbins_Pt_Log, &bins_Pt_Log[0]);
00147     MuonEfficiency_pT_LOGX_sig  = dbe_->book1D("MuonEfficiency_pT_LOGX_sig",  "P^{#mu}_{T} (signal triggered)",  nbins_Pt_Log, &bins_Pt_Log[0]);
00148     MuonEfficiency_pT_LOGX_trig = dbe_->book1D("MuonEfficiency_pT_LOGX_trig", "P^{#mu}_{T} (control triggered)", nbins_Pt_Log, &bins_Pt_Log[0]);
00149 
00150     const int nbins_mass = 200;
00151 
00152     float bins_mass[nbins_mass+1];
00153 
00154     for (int i = 0; i <= nbins_mass; i++) {
00155       double log = logmin + (logmax-logmin)*i/nbins_mass;
00156       bins_mass[i] = std::pow(10.0, log);
00157     }
00158 
00159     DiMuonMassRC       = dbe_->book1D("DiMuonMassRC",      "Invariant Dimuon Mass (Right Charge)", 50, 0., 200.);
00160     DiMuonMassWC       = dbe_->book1D("DiMuonMassWC",      "Invariant Dimuon Mass (Wrong Charge)", 50, 0., 200.);
00161 
00162     DiMuonMassRC_LOGX  = dbe_->book1D("DiMuonMassRC_LOGX", "Invariant Dimuon Mass (Right Charge)", nbins_mass, &bins_mass[0]);
00163     DiMuonMassWC_LOGX  = dbe_->book1D("DiMuonMassWC_LOGX", "Invariant Dimuon Mass (Wrong Charge)", nbins_mass, &bins_mass[0]);
00164 
00165   }
00166 
00167 }
00168 
00169 
00170 void TopHLTDiMuonDQM::beginRun(const edm::Run& r, const edm::EventSetup& context) {
00171 
00172 }
00173 
00174 
00175 void TopHLTDiMuonDQM::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg, const edm::EventSetup& context) {
00176 
00177 }
00178 
00179 
00180 void TopHLTDiMuonDQM::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup ) {
00181 
00182   // ------------------------
00183   //  Global Event Variables
00184   // ------------------------
00185 
00186   vector<string> hltPaths = hltPaths_L3_;
00187 
00188   vector<reco::Particle> Triggered_muons;
00189   reco::MuonCollection    Isolated_muons;
00190   reco::MuonCollection     Matched_muons;
00191 
00192   const int N_TriggerPaths = hltPaths.size();
00193   const int N_SignalPaths  = hltPaths_sig_.size();
00194   const int N_ControlPaths = hltPaths_trig_.size();
00195 
00196   bool Fired_Signal_Trigger[ 10] = {false};
00197   bool Fired_Control_Trigger[10] = {false};
00198 
00199   double DilepMass = 0.;
00200 
00201   double deltaR_Trig  = 1000.;
00202   double deltaR_Reco  =    0.;
00203   double deltaR_Match =    0.;
00204 
00205   int N_iso_mu = 0;
00206 
00207   double vertex_X = 100.;
00208   double vertex_Y = 100.;
00209   double vertex_Z = 100.;
00210 
00211   // -------------------------
00212   //  Analyze Trigger Results
00213   // -------------------------
00214 
00215   edm::Handle<TriggerResults> trigResults;
00216   iEvent.getByLabel(triggerResults_, trigResults);
00217 
00218   if( trigResults.failedToGet() ) {
00219 
00220     //    cout << endl << "-----------------------------" << endl;
00221     //    cout << "--- NO TRIGGER RESULTS !! ---" << endl;
00222     //    cout << "-----------------------------" << endl << endl;
00223 
00224   }
00225 
00226   if( !trigResults.failedToGet() ) {
00227 
00228     const edm::TriggerNames & trigName = iEvent.triggerNames(*trigResults);
00229 
00230     for( unsigned int i_Trig = 0; i_Trig < trigResults->size(); ++i_Trig ) {
00231 
00232       if(trigResults->accept(i_Trig)) {
00233 
00234         // Check for all trigger paths
00235 
00236         for( int i = 0; i < N_TriggerPaths; i++ ) {
00237 
00238           if( trigName.triggerName(i_Trig) == hltPaths[i] )  Trigs->Fill(i);
00239 
00240         }
00241 
00242         // Check for signal & control trigger paths
00243 
00244         for( int j = 0; j < N_SignalPaths; ++j ) {
00245 
00246           if( trigName.triggerName(i_Trig) == hltPaths_sig_[j]  )  Fired_Signal_Trigger[j]  = true;
00247 
00248         }
00249 
00250         for( int k = 0; k < N_ControlPaths; ++k ) {
00251 
00252           if( trigName.triggerName(i_Trig) == hltPaths_trig_[k] )  Fired_Control_Trigger[k] = true;
00253 
00254         }
00255 
00256       }
00257 
00258     }
00259 
00260   }
00261 
00262   // -----------------------
00263   //  Analyze Trigger Event
00264   // -----------------------
00265 
00266   edm::Handle<TriggerEvent> triggerEvent;
00267   iEvent.getByLabel(triggerEvent_, triggerEvent);
00268 
00269   if( triggerEvent.failedToGet() ) {
00270 
00271     //    cout << endl << "---------------------------" << endl;
00272     //    cout << "--- NO TRIGGER EVENT !! ---" << endl;
00273     //    cout << "---------------------------" << endl << endl;
00274 
00275   }
00276 
00277   if( !triggerEvent.failedToGet() ) {
00278 
00279     size_t filterIndex = triggerEvent->filterIndex( triggerFilter_ );
00280     TriggerObjectCollection triggerObjects = triggerEvent->getObjects();
00281 
00282     if( filterIndex >= triggerEvent->sizeFilters() ) {
00283 
00284       //      cout << endl << "------------------------------" << endl;
00285       //      cout << "--- NO FILTERED OBJECTS !! ---" << endl;
00286       //      cout << "------------------------------" << endl << endl;
00287 
00288     }
00289 
00290     if( filterIndex < triggerEvent->sizeFilters() ) {
00291 
00292       const Keys & keys = triggerEvent->filterKeys( filterIndex );
00293 
00294       int N_mu = 0;
00295 
00296       for( size_t j = 0; j < keys.size(); j++ ) {
00297 
00298         TriggerObject foundObject = triggerObjects[keys[j]];
00299 
00300         if(      foundObject.pt()   < muon_pT_cut_  )  continue;
00301         if( abs( foundObject.eta()) > muon_eta_cut_ )  continue;
00302         if( abs( foundObject.particle().pdgId() ) != 13 )  continue;
00303 
00304         ++N_mu;
00305         Triggered_muons.push_back( foundObject.particle() );
00306 
00307       }
00308 
00309       if( Triggered_muons.size() == 2 ) {
00310 
00311         reco::Particle mu1 = Triggered_muons.at(0);
00312         reco::Particle mu2 = Triggered_muons.at(1);
00313 
00314         deltaR_Trig = deltaR( mu1.eta(), mu1.phi(), mu2.eta(), mu2.phi() );
00315         DeltaR_Trig->Fill(deltaR_Trig);
00316 
00317       }
00318 
00319     }
00320 
00321   }
00322 
00323   // ------------------------
00324   //  Analyze Primary Vertex
00325   // ------------------------
00326 
00327   edm::Handle<reco::VertexCollection> vertexs;
00328   iEvent.getByLabel(vertex_, vertexs);
00329 
00330   if( vertexs.failedToGet() ) {
00331 
00332     //      cout << endl << "----------------------------" << endl;
00333     //      cout << "--- NO PRIMARY VERTEX !! ---" << endl;
00334     //      cout << "----------------------------" << endl << endl;
00335 
00336   }
00337 
00338   if( !vertexs.failedToGet() ) {
00339 
00340     reco::Vertex primaryVertex = vertexs->front();
00341 
00342     int numberTracks = primaryVertex.tracksSize();
00343     //      double ndof      = primaryVertex.ndof();
00344     bool fake        = primaryVertex.isFake();
00345 
00346     NTracks->Fill(numberTracks);
00347 
00348     if( !fake && numberTracks > 3 ) {
00349 
00350       vertex_X = primaryVertex.x();
00351       vertex_Y = primaryVertex.y();
00352       vertex_Z = primaryVertex.z();
00353 
00354     }
00355 
00356   }
00357 
00358   // --------------------
00359   //  Analyze RECO Muons
00360   // --------------------
00361 
00362   edm::Handle<reco::MuonCollection> muons;
00363   iEvent.getByLabel(muons_, muons);
00364 
00365   reco::MuonCollection::const_iterator muon;
00366 
00367   if( muons.failedToGet() ) {
00368 
00369     //    cout << endl << "------------------------" << endl;
00370     //    cout << "--- NO RECO MUONS !! ---" << endl;
00371     //    cout << "------------------------" << endl << endl;
00372 
00373   }
00374 
00375   if( !muons.failedToGet() ) {
00376 
00377     NMuons->Fill( muons->size() );
00378 
00379     //    cout << "N_muons : " << muons->size() << endl;
00380 
00381     for(muon = muons->begin(); muon!= muons->end(); ++muon) {
00382 
00383       float N_muons = muons->size();
00384       float Q_muon  = muon->charge();
00385 
00386       NMuons_charge->Fill(N_muons*Q_muon);
00387 
00388       double track_X = 100.;
00389       double track_Y = 100.;
00390       double track_Z = 100.;
00391 
00392       double N_PixelHits   = 0.;
00393       double N_TrackerHits = 0.;
00394 
00395       if( muon->isGlobalMuon() && muon->isTrackerMuon() ) {
00396 
00397         reco::TrackRef track = muon->globalTrack();
00398 
00399         track_X       = track->vx();
00400         track_Y       = track->vy();
00401         track_Z       = track->vz();
00402         N_PixelHits   = track->hitPattern().numberOfValidPixelHits();
00403         N_TrackerHits = track->hitPattern().numberOfValidTrackerHits();
00404 
00405         VxVy_muons->Fill(track_X, track_Y);
00406         Vz_muons->Fill(track_Z);
00407         PixelHits_muons->Fill(N_PixelHits);
00408         TrackerHits_muons->Fill(N_TrackerHits);
00409 
00410       }
00411 
00412       // Vertex and kinematic cuts
00413 
00414       if(          track_X > vertex_X_cut_ )  continue;
00415       if(          track_Y > vertex_Y_cut_ )  continue;
00416       if(          track_Z > vertex_Z_cut_ )  continue;
00417       if(      N_PixelHits <  1.           )  continue;
00418       if(    N_TrackerHits < 11.           )  continue;
00419       if(     muon->pt()   < muon_pT_cut_  )  continue;
00420       if( abs(muon->eta()) > muon_eta_cut_ )  continue;
00421 
00422       reco::MuonIsolation muIso03 = muon->isolationR03();
00423 
00424       double muonCombRelIso = 1.;
00425 
00426       if ( muon->pt() != 0. )
00427         muonCombRelIso = ( muIso03.emEt + muIso03.hadEt + muIso03.hoEt + muIso03.sumPt ) / muon->pt();
00428 
00429       CombRelIso03->Fill( muonCombRelIso );
00430 
00431       if( muonCombRelIso < muon_iso_cut_ ) {
00432 
00433         ++N_iso_mu;
00434         Isolated_muons.push_back(*muon);
00435 
00436       }
00437 
00438     }
00439 
00440     NMuons_iso->Fill(N_iso_mu);
00441 
00442     //    if( Isolated_muons.size() > 1 && Fired_Control_Trigger[0] ) {
00443     if( Isolated_muons.size() > 1 ) {
00444 
00445       // Vertex cut
00446 
00447       if( vertex_X < vertex_X_cut_ && vertex_Y < vertex_Y_cut_ && vertex_Z < vertex_Z_cut_ ) {
00448 
00449         for( int i = 0; i < (static_cast<int>(Isolated_muons.size()) - 1); ++i ) {
00450 
00451           for( int j = i+1; j < static_cast<int>(Isolated_muons.size()); ++j ) {
00452 
00453             reco::MuonCollection::const_reference mu1 = Isolated_muons.at(i);
00454             reco::MuonCollection::const_reference mu2 = Isolated_muons.at(j);
00455 
00456             DilepMass = sqrt( (mu1.energy()+mu2.energy())*(mu1.energy()+mu2.energy())
00457                               - (mu1.px()+mu2.px())*(mu1.px()+mu2.px())
00458                               - (mu1.py()+mu2.py())*(mu1.py()+mu2.py())
00459                               - (mu1.pz()+mu2.pz())*(mu1.pz()+mu2.pz())
00460                               );
00461 
00462             if( DilepMass < 1. ) {
00463 
00464               if( i > 0 ) {
00465 
00466                 Isolated_muons.erase(Isolated_muons.begin()+i);
00467                 --i;
00468 
00469               }
00470 
00471               continue;
00472 
00473             }
00474 
00475             // Opposite muon charges -> Right Charge (RC)
00476 
00477             if( mu1.charge()*mu2.charge() < 0. ) {
00478 
00479               DiMuonMassRC->Fill( DilepMass );
00480               DiMuonMassRC_LOGX->Fill( DilepMass );
00481 
00482               if( DilepMass > MassWindow_down_ && DilepMass < MassWindow_up_ ) {
00483 
00484                 PtMuons->Fill(  mu1.pt()  );
00485                 PtMuons->Fill(  mu2.pt()  );
00486                 PtMuons_LOGX->Fill( mu1.pt() );
00487                 PtMuons_LOGX->Fill( mu2.pt() );
00488                 EtaMuons->Fill( mu1.eta() );
00489                 EtaMuons->Fill( mu2.eta() );
00490                 PhiMuons->Fill( mu1.phi() );
00491                 PhiMuons->Fill( mu2.phi() );
00492 
00493                 DeltaEtaMuonsRC->Fill(mu1.eta()-mu2.eta());
00494                 DeltaPhiMuonsRC->Fill( deltaPhi(mu1.phi(),mu2.phi()) );
00495 
00496                 // Determinating relative trigger efficiencies
00497 
00498                 for( int k = 0; k < N_SignalPaths; ++k ) {
00499 
00500                   if( Fired_Signal_Trigger[k] && Fired_Control_Trigger[k] )  TriggerEfficiencies_sig->Fill(k);
00501 
00502                   if( Fired_Control_Trigger[k] )  TriggerEfficiencies_trig->Fill(k);
00503 
00504                 }
00505 
00506                 // Trigger object matching
00507 
00508                 int N_Match = 0;
00509                 double   DR = 0.1;
00510 
00511                 if( Isolated_muons.size() == 2 && Triggered_muons.size() > 0 ) {
00512 
00513                   deltaR_Reco = deltaR( mu1.eta(), mu1.phi(), mu2.eta(), mu2.phi() );
00514                   DeltaR_Reco->Fill(deltaR_Reco);
00515 
00516                   if( deltaR_Reco > 2.*DR && deltaR_Trig > 2.*DR ) {
00517 
00518                     for( int i = 0; i < static_cast<int>(Isolated_muons.size()); ++i ) {
00519 
00520                       for( int j = 0; j < static_cast<int>(Triggered_muons.size()); ++j ) {
00521 
00522                         reco::Particle & Trigger_mu = Triggered_muons.at(j);
00523                         reco::Muon     &    Reco_mu =  Isolated_muons.at(i);
00524 
00525                         deltaR_Match = deltaR( Trigger_mu.eta(), Trigger_mu.phi(), Reco_mu.eta(), Reco_mu.phi() );
00526                         DeltaR_Match->Fill(deltaR_Match);
00527 
00528                         if( deltaR_Match < DR) {
00529 
00530                           ++N_Match;
00531                           Matched_muons.push_back(Reco_mu);
00532 
00533                         }
00534 
00535                       }
00536 
00537                     }
00538 
00539                     Trigger_Match->Fill(N_Match);
00540 
00541                   }
00542 
00543                 }
00544 
00545 
00546                 // Muon Tag & Probe Efficiency
00547 
00548                 if( Matched_muons.size() == 1 ) {
00549 
00550                   reco::MuonCollection::const_reference matched_mu1 = Matched_muons.at(0);
00551 
00552                   MuonEfficiency_pT_trig->Fill( matched_mu1.pt() );
00553                   MuonEfficiency_pT_LOGX_trig->Fill( matched_mu1.pt() );
00554                   MuonEfficiency_eta_trig->Fill(matched_mu1.eta());
00555                   MuonEfficiency_phi_trig->Fill(matched_mu1.phi());
00556 
00557                 }
00558 
00559                 if( Matched_muons.size() == 2 ) {
00560 
00561                   reco::MuonCollection::const_reference matched_mu1 = Matched_muons.at(0);
00562                   reco::MuonCollection::const_reference matched_mu2 = Matched_muons.at(1);
00563 
00564                   MuonEfficiency_pT_trig->Fill( matched_mu1.pt() );
00565                   MuonEfficiency_pT_trig->Fill( matched_mu2.pt() );
00566                   MuonEfficiency_pT_LOGX_trig->Fill( matched_mu1.pt() );
00567                   MuonEfficiency_pT_LOGX_trig->Fill( matched_mu2.pt() );
00568                   MuonEfficiency_eta_trig->Fill(matched_mu1.eta());
00569                   MuonEfficiency_eta_trig->Fill(matched_mu2.eta());
00570                   MuonEfficiency_phi_trig->Fill(matched_mu1.phi());
00571                   MuonEfficiency_phi_trig->Fill(matched_mu2.phi());
00572 
00573                   MuonEfficiency_pT_sig->Fill( matched_mu1.pt() );
00574                   MuonEfficiency_pT_sig->Fill( matched_mu2.pt() );
00575                   MuonEfficiency_pT_LOGX_sig->Fill( matched_mu1.pt() );
00576                   MuonEfficiency_pT_LOGX_sig->Fill( matched_mu2.pt() );
00577                   MuonEfficiency_eta_sig->Fill(matched_mu1.eta());
00578                   MuonEfficiency_eta_sig->Fill(matched_mu2.eta());
00579                   MuonEfficiency_phi_sig->Fill(matched_mu1.phi());
00580                   MuonEfficiency_phi_sig->Fill(matched_mu2.phi());
00581 
00582                 }
00583 
00584               }
00585 
00586             }
00587 
00588             // Same muon charges -> Wrong Charge (WC)
00589 
00590             if( mu1.charge()*mu2.charge() > 0. ) {
00591 
00592               DiMuonMassWC->Fill( DilepMass );
00593               DiMuonMassWC_LOGX->Fill( DilepMass );
00594 
00595               if( DilepMass > MassWindow_down_ && DilepMass < MassWindow_up_ ) {
00596 
00597                 DeltaEtaMuonsWC->Fill( mu1.eta()-mu2.eta() );
00598                 DeltaPhiMuonsWC->Fill( deltaPhi(mu1.phi(),mu2.phi()) );
00599 
00600               }
00601 
00602             }
00603 
00604           }
00605 
00606         }
00607 
00608       }
00609 
00610     }
00611 
00612   }
00613 
00614 }
00615 
00616 
00617 void TopHLTDiMuonDQM::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg, const edm::EventSetup& context) {
00618 
00619 }
00620 
00621 
00622 void TopHLTDiMuonDQM::endRun(const edm::Run& r, const edm::EventSetup& context) {
00623 
00624 }
00625 
00626 
00627 void TopHLTDiMuonDQM::endJob() {
00628 
00629 }