CMS 3D CMS Logo

Public Member Functions | Protected Attributes

EfficiencyAnalyzer Class Reference

#include <EfficiencyAnalyzer.h>

Inheritance diagram for EfficiencyAnalyzer:
MuonAnalyzerBase

List of all members.

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
void beginJob (DQMStore *dbe)
 Inizialize parameters for histo binning.
 EfficiencyAnalyzer (const edm::ParameterSet &pset, MuonServiceProxy *theService)
virtual ~EfficiencyAnalyzer ()

Protected Attributes

int etaBin_
double etaMax_
double etaMin_
MonitorElementh_allProbes_barrel_pt
MonitorElementh_allProbes_endcap_pt
MonitorElementh_allProbes_eta
MonitorElementh_allProbes_hp_eta
MonitorElementh_allProbes_phi
MonitorElementh_allProbes_pt
MonitorElementh_failProbes_TightMu_eta
MonitorElementh_failProbes_TightMu_phi
MonitorElementh_failProbes_TightMu_pt
MonitorElementh_passProbes_TightMu_barrel_pt
MonitorElementh_passProbes_TightMu_endcap_pt
MonitorElementh_passProbes_TightMu_eta
MonitorElementh_passProbes_TightMu_hp_eta
MonitorElementh_passProbes_TightMu_phi
MonitorElementh_passProbes_TightMu_pt
std::string metname
edm::ParameterSet parameters
int phiBin_
double phiMax_
double phiMin_
int ptBin_
double ptMax_
double ptMin_
MonitorElementtest_TightMu_Minv
edm::InputTag theMuonCollectionLabel
edm::InputTag theSTACollectionLabel
edm::InputTag theTrackCollectionLabel

Detailed Description

Class EfficiencyAnalyzer

DQM monitoring for dimuon mass

Author: S.Folgueras, A. Calderon

Definition at line 26 of file EfficiencyAnalyzer.h.


Constructor & Destructor Documentation

EfficiencyAnalyzer::EfficiencyAnalyzer ( const edm::ParameterSet pset,
MuonServiceProxy *  theService 
)

Definition at line 38 of file EfficiencyAnalyzer.cc.

References parameters.

EfficiencyAnalyzer::~EfficiencyAnalyzer ( ) [virtual]

Definition at line 42 of file EfficiencyAnalyzer.cc.

{ }

Member Function Documentation

void EfficiencyAnalyzer::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)

to be read from output as "generalTracks"

Definition at line 98 of file EfficiencyAnalyzer.cc.

References SiPixelRawToDigiRegional_cfi::beamSpot, reco::MuonIsolation::emEt, MonitorElement::Fill(), edm::Event::getByLabel(), h_allProbes_barrel_pt, h_allProbes_endcap_pt, h_allProbes_eta, h_allProbes_hp_eta, h_allProbes_phi, h_allProbes_pt, h_passProbes_TightMu_barrel_pt, h_passProbes_TightMu_endcap_pt, h_passProbes_TightMu_eta, h_passProbes_TightMu_hp_eta, h_passProbes_TightMu_phi, h_passProbes_TightMu_pt, reco::MuonIsolation::hadEt, edm::HandleBase::isValid(), LogTrace, metname, patZpeak::muons, reco::BeamSpot::position(), reco::MuonIsolation::sumPt, test_TightMu_Minv, theMuonCollectionLabel, theTrackCollectionLabel, and testEve_cfg::tracks.

                                                                                    {

  LogTrace(metname)<<"[EfficiencyAnalyzer] Analyze the mu in different eta regions";
  
  edm::Handle<reco::MuonCollection> muons;
  iEvent.getByLabel(theMuonCollectionLabel, muons);


  edm::Handle<reco::TrackCollection> tracks;
  iEvent.getByLabel(theTrackCollectionLabel,tracks); 


  reco::BeamSpot beamSpot;
  Handle<reco::BeamSpot> beamSpotHandle;
  iEvent.getByLabel("offlineBeamSpot", beamSpotHandle);
  beamSpot = *beamSpotHandle;

  if(!muons.isValid()) return;


  // Loop on muon collection
  TLorentzVector Mu1, Mu2;

  bool isMB = false;
  bool isME = false;

  for (reco::MuonCollection::const_iterator recoMu1 = muons->begin(); recoMu1!=muons->end(); ++recoMu1) {

    LogTrace(metname)<<"[EfficiencyAnalyzer] loop over first muons" << endl;

    //--- Define combined isolation
    reco::MuonIsolation Iso_muon = recoMu1->isolationR03();
    float combIso = (Iso_muon.emEt + Iso_muon.hadEt + Iso_muon.sumPt);  
    
    //--- Is Global Muon 
    if (!recoMu1->isGlobalMuon()) continue;

      // get the track combinig the information from both the Tracker and the Spectrometer
      reco::TrackRef recoCombinedGlbTrack1 = recoMu1->combinedMuon();    
      float muPt1 = recoCombinedGlbTrack1->pt();
      Mu1.SetPxPyPzE(recoCombinedGlbTrack1->px(), recoCombinedGlbTrack1->py(),recoCombinedGlbTrack1->pz(), recoCombinedGlbTrack1->p());
      

      //--- Define if it is a tight muon
      if (recoMu1->isGlobalMuon() && recoMu1->isTrackerMuon() && recoMu1->combinedMuon()->normalizedChi2()<10. 
          && recoMu1->combinedMuon()->hitPattern().numberOfValidMuonHits()>0 && fabs(recoMu1->combinedMuon()->dxy(beamSpot.position()))<0.2 && recoMu1->combinedMuon()->hitPattern().numberOfValidPixelHits()>0 && recoMu1->numberOfMatches() > 1) {

        //-- is isolated muon
        if (muPt1 > 15  && (combIso/muPt1) < 0.1 ) {


          for (reco::MuonCollection::const_iterator recoMu2 = muons->begin(); recoMu2!=muons->end(); ++recoMu2){ 
          
            LogTrace(metname)<<"[EfficiencyAnalyzer] loop over second muon" <<endl;
          
            if (recoMu2 == recoMu1) continue;

            
            if (recoMu2->eta() < 1.479 )  isMB = true;
            if (recoMu2->eta() >= 1.479 ) isME = true;


            //--> should we apply track quality cuts??? 
            Mu2.SetPxPyPzE(recoMu2->px(), recoMu2->py(), recoMu2->pz(), recoMu2->p());

            if (!recoMu2->isTrackerMuon()) continue;

            if ( recoMu2->pt() < 5 ) continue;

            if ( (recoMu1->charge())*(recoMu2->charge()) > 0 ) continue; 

            float Minv = (Mu1+Mu2).M();
            
            if ( Minv < 70 ||  Minv > 110 ) continue; 

            h_allProbes_pt->Fill(recoMu2->pt());
            h_allProbes_eta->Fill(recoMu2->eta());
            h_allProbes_phi->Fill(recoMu2->phi());


            if (isMB) h_allProbes_barrel_pt->Fill(recoMu2->pt());
            if (isME) h_allProbes_endcap_pt->Fill(recoMu2->pt());
            if(recoMu2->pt() > 20 ) h_allProbes_hp_eta->Fill(recoMu2->eta());


            test_TightMu_Minv->Fill(Minv);

         
            // Probes passing the tight muon criteria 

            if (recoMu2->isGlobalMuon() && recoMu2->isTrackerMuon() && recoMu2->combinedMuon()->normalizedChi2()<10. && recoMu2->combinedMuon()->hitPattern().numberOfValidMuonHits()>0 && fabs(recoMu2->combinedMuon()->dxy(beamSpot.position()))<0.2 && recoMu2->combinedMuon()->hitPattern().numberOfValidPixelHits()>0 && recoMu2->numberOfMatches() > 1) { 
                 
                 h_passProbes_TightMu_pt->Fill(recoMu2->pt());
                 h_passProbes_TightMu_eta->Fill(recoMu2->eta());
                 h_passProbes_TightMu_phi->Fill(recoMu2->phi());

                 if (isMB) h_passProbes_TightMu_barrel_pt->Fill(recoMu2->pt());
                 if (isME) h_passProbes_TightMu_endcap_pt->Fill(recoMu2->pt());
                 if( recoMu2->pt() > 20 ) h_passProbes_TightMu_hp_eta->Fill(recoMu2->eta());
       
            }
          }
              
        }
      }
  }
}
void EfficiencyAnalyzer::beginJob ( DQMStore dbe) [virtual]

Inizialize parameters for histo binning.

Implements MuonAnalyzerBase.

Definition at line 44 of file EfficiencyAnalyzer.cc.

References DQMStore::book1D(), gather_cfg::cout, etaBin_, etaMax_, etaMin_, edm::ParameterSet::getParameter(), h_allProbes_barrel_pt, h_allProbes_endcap_pt, h_allProbes_eta, h_allProbes_hp_eta, h_allProbes_phi, h_allProbes_pt, h_passProbes_TightMu_barrel_pt, h_passProbes_TightMu_endcap_pt, h_passProbes_TightMu_eta, h_passProbes_TightMu_hp_eta, h_passProbes_TightMu_phi, h_passProbes_TightMu_pt, LogTrace, metname, parameters, phiBin_, phiMax_, phiMin_, ptBin_, ptMax_, ptMin_, DQMStore::setCurrentFolder(), test_TightMu_Minv, theMuonCollectionLabel, and theTrackCollectionLabel.

                                                {
#ifdef DEBUG
  cout << "[EfficiencyAnalyzer] Parameters initialization" <<endl;
#endif
  metname = "EfficiencyAnalyzer";
  LogTrace(metname)<<"[EfficiencyAnalyzer] Parameters initialization";
  dbe->setCurrentFolder("Muons/EfficiencyAnalyzer");  
  
  theMuonCollectionLabel = parameters.getParameter<edm::InputTag>("MuonCollection");
  theTrackCollectionLabel = parameters.getParameter<edm::InputTag>("TrackCollection");


  ptBin_ = parameters.getParameter<int>("ptBin");
  ptMin_ = parameters.getParameter<double>("ptMin");
  ptMax_ = parameters.getParameter<double>("ptMax");

  etaBin_ = parameters.getParameter<int>("etaBin");
  etaMin_ = parameters.getParameter<double>("etaMin");
  etaMax_ = parameters.getParameter<double>("etaMax");

  phiBin_ = parameters.getParameter<int>("phiBin");
  phiMin_ = parameters.getParameter<double>("phiMin");
  phiMax_ = parameters.getParameter<double>("phiMax");



  test_TightMu_Minv  = dbe->book1D("test_TightMu_Minv"  ,"Minv",50,70,110);

  h_allProbes_pt = dbe->book1D("allProbes_pt","All Probes Pt", ptBin_, ptMin_, ptMax_);
  h_allProbes_barrel_pt = dbe->book1D("allProbes_barrel_pt","Barrel: all Probes Pt", ptBin_, ptMin_, ptMax_);
  h_allProbes_endcap_pt = dbe->book1D("allProbes_endcap_pt","Endcap: all Probes Pt", ptBin_, ptMin_, ptMax_);
  h_allProbes_eta = dbe->book1D("allProbes_eta","All Probes Eta", etaBin_, etaMin_, etaMax_);
  h_allProbes_hp_eta = dbe->book1D("allProbes_hp_eta","High Pt all Probes Eta", etaBin_, etaMin_, etaMax_);
  h_allProbes_phi = dbe->book1D("allProbes_phi","All Probes Phi", phiBin_, phiMin_, phiMax_);

  h_passProbes_TightMu_pt = dbe->book1D("passProbes_TightMu_pt","TightMu Passing Probes Pt", ptBin_ , ptMin_ , ptMax_ );
 h_passProbes_TightMu_barrel_pt = dbe->book1D("passProbes_TightMu_barrel_pt","Barrel: TightMu Passing Probes Pt", ptBin_ , ptMin_ , ptMax_ );
 h_passProbes_TightMu_endcap_pt = dbe->book1D("passProbes_TightMu_endcap_pt","Endcap: TightMu Passing Probes Pt", ptBin_ , ptMin_ , ptMax_ );
  h_passProbes_TightMu_eta = dbe->book1D("passProbes_TightMu_eta","TightMu Passing Probes #eta", etaBin_, etaMin_, etaMax_);
  h_passProbes_TightMu_hp_eta = dbe->book1D("passProbes_TightMu_hp_eta","High Pt TightMu Passing Probes #eta", etaBin_, etaMin_, etaMax_);
  h_passProbes_TightMu_phi = dbe->book1D("passProbes_TightMu_phi","TightMu Passing Probes #phi", phiBin_, phiMin_, phiMax_);



  /*h_failProbes_TightMu_pt = dbe->book1D("failProbes_TightMu_pt","TightMu Failling Probes Pt", ptBin_ , ptMin_ , ptMax_ );
  h_failProbes_TightMu_eta = dbe->book1D("failProbes_TightMu_eta","TightMu Failling Probes #eta", etaBin_, etaMin_, etaMax_);
  h_failProbes_TightMu_phi = dbe->book1D("failProbes_TightMu_phi","TightMu Failling Probes #phi", phiBin_, phiMin_, phiMax_);
  */

#ifdef DEBUG
  cout << "[EfficiencyAnalyzer] Parameters initialization DONE" <<endl;
#endif
}

Member Data Documentation

int EfficiencyAnalyzer::etaBin_ [protected]

Definition at line 53 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

double EfficiencyAnalyzer::etaMax_ [protected]

Definition at line 61 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

double EfficiencyAnalyzer::etaMin_ [protected]

Definition at line 60 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

Definition at line 78 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 79 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 80 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 81 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 82 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 77 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 74 of file EfficiencyAnalyzer.h.

Definition at line 75 of file EfficiencyAnalyzer.h.

Definition at line 73 of file EfficiencyAnalyzer.h.

Definition at line 67 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 68 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 69 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 70 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 71 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 66 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

std::string EfficiencyAnalyzer::metname [protected]

Definition at line 44 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 41 of file EfficiencyAnalyzer.h.

Referenced by beginJob(), and EfficiencyAnalyzer().

int EfficiencyAnalyzer::phiBin_ [protected]

Definition at line 54 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

double EfficiencyAnalyzer::phiMax_ [protected]

Definition at line 64 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

double EfficiencyAnalyzer::phiMin_ [protected]

Definition at line 63 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

int EfficiencyAnalyzer::ptBin_ [protected]

Definition at line 55 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

double EfficiencyAnalyzer::ptMax_ [protected]

Definition at line 58 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

double EfficiencyAnalyzer::ptMin_ [protected]

Definition at line 57 of file EfficiencyAnalyzer.h.

Referenced by beginJob().

Definition at line 85 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 48 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().

Definition at line 47 of file EfficiencyAnalyzer.h.

Definition at line 49 of file EfficiencyAnalyzer.h.

Referenced by analyze(), and beginJob().