CMS 3D CMS Logo

Public Types | Public Member Functions | Private Member Functions | Private Attributes

MuonPFAnalyzer Class Reference

#include <MuonPFAnalyzer.h>

Inheritance diagram for MuonPFAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

List of all members.

Public Types

typedef std::vector< RecoGenPairRecoGenCollection
typedef std::pair< const
reco::Muon *, const
reco::GenParticle * > 
RecoGenPair

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 Perform the PF - TUNEP muon analysis.
virtual void beginRun (edm::Run const &, edm::EventSetup const &)
 Initialize an book plots.
 MuonPFAnalyzer (const edm::ParameterSet &)
 Constructor.
 ~MuonPFAnalyzer ()
 Destructor.

Private Member Functions

void bookHistos (const std::string &group)
float combRelIso (const reco::Muon *muon)
float fDeltaPhi (float phi1, float phi2)
void fillInRange (MonitorElement *plot, int nAxis, double x, double y=0)
MonitorElementgetPlot (const std::string &group, const std::string &type)
const reco::Vertex getPrimaryVertex (edm::Handle< reco::VertexCollection > &vertex, edm::Handle< reco::BeamSpot > &beamSpot)
int muonTrackType (const reco::Muon *muon, bool usePF)
void recoToGenMatch (edm::Handle< reco::MuonCollection > &reco, edm::Handle< reco::GenParticleCollection > &gen)
void setCodeLabels (MonitorElement *plot, int nAxis)

Private Attributes

edm::InputTag theBeamSpotLabel
DQMStoretheDbe
std::string theFolder
edm::InputTag theGenLabel
double theHighPtTh
double theIsoCut
std::vector< std::string > theMuonKinds
std::map< std::string,
std::map< std::string,
MonitorElement * > > 
thePlots
RecoGenCollection theRecoGen
double theRecoGenR
edm::InputTag theRecoLabel
bool theRunOnMC
edm::InputTag theVertexLabel

Detailed Description

Definition at line 36 of file MuonPFAnalyzer.h.


Member Typedef Documentation

Definition at line 41 of file MuonPFAnalyzer.h.

typedef std::pair<const reco::Muon*, const reco::GenParticle*> MuonPFAnalyzer::RecoGenPair

Definition at line 40 of file MuonPFAnalyzer.h.


Constructor & Destructor Documentation

MuonPFAnalyzer::MuonPFAnalyzer ( const edm::ParameterSet pSet) [explicit]

Constructor.

Definition at line 45 of file MuonPFAnalyzer.cc.

References edm::ParameterSet::getParameter(), and LogTrace.

{

  LogTrace("MuonPFAnalyzer") << 
    "[MuonPFAnalyzer] Initializing configuration from parameterset.\n";

  theGenLabel      = pSet.getParameter<InputTag>("inputTagGenParticles");  
  theRecoLabel     = pSet.getParameter<InputTag>("inputTagMuonReco");  
  theBeamSpotLabel = pSet.getParameter<InputTag>("inputTagBeamSpot");  
  theVertexLabel   = pSet.getParameter<InputTag>("inputTagVertex");  
  
  theHighPtTh   = pSet.getParameter<double>("highPtThreshold");
  theRecoGenR   = pSet.getParameter<double>("recoGenDeltaR");
  theIsoCut     = pSet.getParameter<double>("relCombIsoCut");
  theRunOnMC    = pSet.getParameter<bool>("runOnMC");

  theFolder = pSet.getParameter<string>("folder");
  
  theMuonKinds.push_back("");          // all TUNEP/PF muons
  theMuonKinds.push_back("Tight");     // tight TUNEP/PF muons
  theMuonKinds.push_back("TightIso");  // tight/iso TUNEP/PF muons
 

}
MuonPFAnalyzer::~MuonPFAnalyzer ( )

Destructor.

Definition at line 72 of file MuonPFAnalyzer.cc.

References LogTrace.

{

  LogTrace("MuonPFAnalyzer") << 
    "[MuonPFAnalyzer] Destructor called.\n";
  
}

Member Function Documentation

void MuonPFAnalyzer::analyze ( const edm::Event event,
const edm::EventSetup context 
) [virtual]

Perform the PF - TUNEP muon analysis.

Implements edm::EDAnalyzer.

Definition at line 111 of file MuonPFAnalyzer.cc.

References SiPixelRawToDigiRegional_cfi::beamSpot, eta(), dataDML::fillInRange(), reco::Muon::innerTrack(), edm::Ref< C, T, F >::isNull(), reco::Muon::isPFMuon(), muon::isTightMuon(), metsig::muon, patZpeak::muons, reco::LeafCandidate::p4(), phi, reco::LeafCandidate::pt(), and reco::Muon::tunePMuonBestTrack().

{
  
  Handle<reco::MuonCollection> muons;
  event.getByLabel(theRecoLabel, muons);

  Handle<GenParticleCollection> genMuons;
  event.getByLabel(theGenLabel, genMuons);

  Handle<BeamSpot> beamSpot;
  event.getByLabel(theBeamSpotLabel, beamSpot);

  Handle<VertexCollection> vertex;
  event.getByLabel(theVertexLabel, vertex);

  const Vertex primaryVertex = getPrimaryVertex(vertex, beamSpot);

  recoToGenMatch(muons, genMuons);
 
  RecoGenCollection::const_iterator recoGenIt  = theRecoGen.begin();
  RecoGenCollection::const_iterator recoGenEnd = theRecoGen.end();
    
  for (;recoGenIt!=recoGenEnd;++recoGenIt) 
    {
    
      const Muon *muon = recoGenIt->first;
      TrackRef tunePTrack = muon->tunePMuonBestTrack();

     const GenParticle *genMuon = recoGenIt->second;
          
      vector<string>::const_iterator kindIt  = theMuonKinds.begin();
      vector<string>::const_iterator kindEnd = theMuonKinds.end();

      for (;kindIt!=kindEnd;++kindIt) 
        { 

          const string & kind = (*kindIt);

          if (kind.find("Tight") != string::npos && 
              !muon::isTightMuon((*muon), primaryVertex)) continue;
          
          if (kind.find("Iso") != string::npos &&
              combRelIso(muon) > theIsoCut) continue;
          
          if (theRunOnMC && genMuon && !muon->innerTrack().isNull() ) // has matched gen muon
            {

              if (!tunePTrack.isNull())
                { 

                  string group = "TUNEP" + kind;

                  float pt  = tunePTrack->pt();
                  float phi = tunePTrack->phi();
                  float eta = tunePTrack->eta();

                  float genPt  = genMuon->pt();
                  float genPhi = genMuon->p4().phi();
                  float genEta = genMuon->p4().eta();

                  float dPtOverPt = (pt / genPt) - 1;
                  
                  if (pt < theHighPtTh)
                    {

                      fillInRange(getPlot(group,"code"),1,muonTrackType(muon, false));
                      fillInRange(getPlot(group,"deltaPtOverPt"),1,dPtOverPt);
                    }
                  else
                    {
                      fillInRange(getPlot(group,"codeHighPt"),1,muonTrackType(muon, false));
                      fillInRange(getPlot(group,"deltaPtOverPtHighPt"),1,dPtOverPt);
                    }
    
                  fillInRange(getPlot(group,"deltaPt"),1,(pt - genPt));
                  fillInRange(getPlot(group,"deltaPhi"),1,fDeltaPhi(genPhi,phi));
                  fillInRange(getPlot(group,"deltaEta"),1,genEta - eta);
          
                }
              
              if (muon->isPFMuon()) 
                {
                  
                  string group = "PF" + kind;
                  
                  // Assumes that default in muon is PF
                  float pt  = muon->pt();
                  float phi = muon->p4().phi();
                  float eta = muon->p4().eta();
                  
                  float genPt  = genMuon->pt();
                  float genPhi = genMuon->p4().phi();
                  float genEta = genMuon->p4().eta();

                  float dPtOverPt = (pt / genPt) - 1;
                  
                  if (pt < theHighPtTh)
                    {
                      fillInRange(getPlot(group,"code"),1,muonTrackType(muon, true));
                      fillInRange(getPlot(group,"deltaPtOverPt"),1,dPtOverPt);
                    }
                  else
                    { 
                      fillInRange(getPlot(group,"codeHighPt"),1,muonTrackType(muon, true));
                      fillInRange(getPlot(group,"deltaPtOverPtHighPt"),1,dPtOverPt);
                    }
                  
                  
                  fillInRange(getPlot(group,"deltaPt"),1,pt - genPt);
                  fillInRange(getPlot(group,"deltaPhi"),1,fDeltaPhi(genPhi,phi));
                  fillInRange(getPlot(group,"deltaEta"),1,genEta - eta);
                  
                }
            
            }

        

            if (muon->isPFMuon() && !tunePTrack.isNull() &&          
                !muon->innerTrack().isNull()) // Compare PF with TuneP + Tracker 
              {                               // No gen matching needed
        
              string group = "PFvsTUNEP" + kind;

              float pt  = tunePTrack->pt();
              float phi = tunePTrack->phi();
              float eta = tunePTrack->eta();
              
              // Assumes that default in muon is PF
              float pfPt  = muon->pt();
              float pfPhi = muon->p4().phi();
              float pfEta = muon->p4().eta();
              float dPtOverPt = (pfPt / pt) - 1; // TUNEP vs PF pt used as denum.
             

              if (pt < theHighPtTh) 
                {
                  fillInRange(getPlot(group,"code"),2,
                              muonTrackType(muon, false),muonTrackType(muon, true));
                  fillInRange(getPlot(group,"deltaPtOverPt"),1,dPtOverPt);
                }
              else 
                {
                  fillInRange(getPlot(group,"codeHighPt"),
                              2,muonTrackType(muon, false),muonTrackType(muon, true));
                  fillInRange(getPlot(group,"deltaPtOverPtHighPt"),1,dPtOverPt);
                }
              
              fillInRange(getPlot(group,"deltaPt"),1,pfPt - pt);
              fillInRange(getPlot(group,"deltaPhi"),1,fDeltaPhi(pfPhi,phi));
              fillInRange(getPlot(group,"deltaEta"),1,pfEta - eta);
              

              if (theRunOnMC && genMuon) // has a matched gen muon
                
                {
                  
                  float genPt     = genMuon->pt();
                  float dPtOverPtGen = (pt / genPt) - 1;
                  float dPtOverPtGenPF = (pfPt / genPt) - 1;
                  
                  if (pt < theHighPtTh) 
                    {
                      fillInRange(getPlot(group,"deltaPtOverPtPFvsTUNEP"),
                                  2,dPtOverPtGen,dPtOverPtGenPF);
                    }
                  else 
                    {
                      fillInRange(getPlot(group,"deltaPtOverPtHighPtPFvsTUNEP"),
                                  2,dPtOverPtGen,dPtOverPtGenPF);
                    }
                }                 
              
              }
            
        }
      
    }

}
void MuonPFAnalyzer::beginRun ( edm::Run const &  ,
edm::EventSetup const &   
) [virtual]

Initialize an book plots.

Reimplemented from edm::EDAnalyzer.

Definition at line 84 of file MuonPFAnalyzer.cc.

References bookHistos(), LogTrace, and cppFunctionSkipper::operator.

                                                                   {

  LogTrace("MuonPFAnalyzer") << 
    "[MuonPFAnalyzer] Booking histograms.\n";

  //Set up DAQ
  theDbe = 0;
  theDbe = edm::Service<DQMStore>().operator->();
  theDbe->cd();

  if(theRunOnMC)
    {
      bookHistos("PF");
      bookHistos("PFTight");
      bookHistos("PFTightIso");
      bookHistos("TUNEP");
      bookHistos("TUNEPTight");
      bookHistos("TUNEPTightIso");
    }
 
    bookHistos("PFvsTUNEP");
    bookHistos("PFvsTUNEPTight");
    bookHistos("PFvsTUNEPTightIso");
  

}
void MuonPFAnalyzer::bookHistos ( const std::string &  group) [private]
float MuonPFAnalyzer::combRelIso ( const reco::Muon muon) [inline, private]
float MuonPFAnalyzer::fDeltaPhi ( float  phi1,
float  phi2 
) [inline, private]

Definition at line 405 of file MuonPFAnalyzer.cc.

References funct::cos().

                                                             {
  
  float fPhiDiff = fabs(acos(cos(phi1-phi2)));
  return fPhiDiff;
  
}
void MuonPFAnalyzer::fillInRange ( MonitorElement plot,
int  nAxis,
double  x,
double  y = 0 
) [private]

Definition at line 439 of file MuonPFAnalyzer.cc.

References MonitorElement::Fill(), MonitorElement::getTH1(), timingPdfMaker::histo, i, edm::max(), edm::min(), relativeConstraints::value, x, and detailsBasic3DVector::y.

{

  TH1 * histo =  plot->getTH1();
  
  TAxis *axis[2] = {0, 0};
  axis[0] = histo->GetXaxis();
  if (nAxis == 2)
    axis[1] = histo->GetYaxis();

  double value[2] = {0, 0};
  value[0] = x;
  value[1] = y;

  for (int i = 0;i<nAxis;++i)
    {
      double min = axis[i]->GetXmin();
      double max = axis[i]->GetXmax();

      if (value[i] <= min)
        value[i] = axis[i]->GetBinCenter(1);

      if (value[i] >= max)
        value[i] = axis[i]->GetBinCenter(axis[i]->GetNbins());
    }

  if (nAxis == 2)
    plot->Fill(value[0],value[1]);
  else
    plot->Fill(value[0]);
  
}
MonitorElement * MuonPFAnalyzer::getPlot ( const std::string &  group,
const std::string &  type 
) [private]

Definition at line 371 of file MuonPFAnalyzer.cc.

References LogTrace.

                                                              {

  map<string,map<string,MonitorElement *> >::iterator groupIt = thePlots.find(group);
  if (groupIt == thePlots.end()) {
    LogTrace("MuonPFAnalyzer") << "[MuonPFAnalyzer] GROUP : " << group 
                               << " is not a valid plot group. Returning 0.\n";
    return 0;
  }
  
  map<string,MonitorElement *>::iterator typeIt = groupIt->second.find(type);
  if (typeIt == groupIt->second.end()) {
    LogTrace("MuonPFAnalyzer") << "[MuonPFAnalyzer] TYPE : " << type 
                               << " is not a valid type for GROUP : " << group 
                               << ". Returning 0.\n";
    return 0;
  }
  
  return typeIt->second;

} 
const reco::Vertex MuonPFAnalyzer::getPrimaryVertex ( edm::Handle< reco::VertexCollection > &  vertex,
edm::Handle< reco::BeamSpot > &  beamSpot 
) [private]

Definition at line 554 of file MuonPFAnalyzer.cc.

References edm::HandleBase::isValid().

{

  Vertex::Point posVtx;
  Vertex::Error errVtx;

  bool hasPrimaryVertex = false;

  if (vertex.isValid())
    {

      vector<Vertex>::const_iterator vertexIt  = vertex->begin();
      vector<Vertex>::const_iterator vertexEnd = vertex->end();

      for (;vertexIt!=vertexEnd;++vertexIt) 
        {
          if (vertexIt->isValid() && 
              !vertexIt->isFake()) 
            {
              posVtx = vertexIt->position();
              errVtx = vertexIt->error();
              hasPrimaryVertex = true;        
              break;
            }
        }
    }

  if ( !hasPrimaryVertex ) {

    LogInfo("MuonPFAnalyzer") << 
      "[MuonPFAnalyzer] PrimaryVertex not found, use BeamSpot position instead.\n";

    posVtx = beamSpot->position();
    errVtx(0,0) = beamSpot->BeamWidthX();
    errVtx(1,1) = beamSpot->BeamWidthY();
    errVtx(2,2) = beamSpot->sigmaZ();
    
  }

  const Vertex primaryVertex(posVtx,errVtx);

  return primaryVertex;

}
int MuonPFAnalyzer::muonTrackType ( const reco::Muon muon,
bool  usePF 
) [private]

Definition at line 473 of file MuonPFAnalyzer.cc.

References reco::Muon::muonBestTrackType(), and reco::Muon::tunePMuonBestTrackType().

                                                               {

  switch ( usePF ? muon->muonBestTrackType() : muon->tunePMuonBestTrackType() ) {
  case Muon::InnerTrack :
    return 0;
  case Muon::OuterTrack :
    return 1;
  case Muon::CombinedTrack :
    return 2;
  case Muon::TPFMS :
    return 3;
  case Muon::Picky :
    return 4;
  case Muon::DYT :
    return 5;
  case Muon::None :
    return 6;
  }

  return 6;

}
void MuonPFAnalyzer::recoToGenMatch ( edm::Handle< reco::MuonCollection > &  reco,
edm::Handle< reco::GenParticleCollection > &  gen 
) [private]

Definition at line 497 of file MuonPFAnalyzer.cc.

References abs, deltaR(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, and edm::HandleBase::isValid().

{

  theRecoGen.clear();  

  if (muons.isValid())
    {

      MuonCollection::const_iterator muonIt  = muons->begin();
      MuonCollection::const_iterator muonEnd = muons->end();

      for(; muonIt!=muonEnd; ++muonIt) 
        {
      
          float bestDR = 999.;
          const GenParticle *bestGen = 0;

          if (theRunOnMC && gens.isValid()) 
            {
          
              GenParticleCollection::const_iterator genIt  = gens->begin();
              GenParticleCollection::const_iterator genEnd = gens->end();
      
              for(; genIt!=genEnd; ++genIt) 
                {
          
                  if (abs(genIt->pdgId()) == 13 ) 
                    {
                  
                      float muonPhi = muonIt->phi();
                      float muonEta = muonIt->eta();
                  
                      float genPhi = genIt->phi();
                      float genEta = genIt->eta();
                  
                      float dR = deltaR(muonEta,muonPhi,
                                        genEta,genPhi);
                  
                      if (dR < theRecoGenR && dR < bestDR) 
                        {
                          bestDR = dR;
                          bestGen = &(*genIt);
                        }
                      
                    }   
                  
                }
            }
      
          theRecoGen.push_back(RecoGenPair(&(*muonIt), bestGen));

        }
    }
  
}
void MuonPFAnalyzer::setCodeLabels ( MonitorElement plot,
int  nAxis 
) [private]

Definition at line 413 of file MuonPFAnalyzer.cc.

References MonitorElement::getTH1(), and timingPdfMaker::histo.

{

  TAxis *axis = 0;
  
  TH1 * histo = plot->getTH1();
  if(!histo) return;
  
  if (nAxis==1) 
    axis =histo->GetXaxis();
  else if (nAxis == 2)
    axis =histo->GetYaxis();

  if(!axis) return;

  axis->SetBinLabel(1,"Inner Track");
  axis->SetBinLabel(2,"Outer Track");
  axis->SetBinLabel(3,"Combined");
  axis->SetBinLabel(4,"TPFMS");
  axis->SetBinLabel(5,"Picky");
  axis->SetBinLabel(6,"DYT");
  axis->SetBinLabel(7,"None");

}

Member Data Documentation

Definition at line 90 of file MuonPFAnalyzer.h.

Definition at line 94 of file MuonPFAnalyzer.h.

std::string MuonPFAnalyzer::theFolder [private]

Definition at line 105 of file MuonPFAnalyzer.h.

Definition at line 87 of file MuonPFAnalyzer.h.

double MuonPFAnalyzer::theHighPtTh [private]

Definition at line 99 of file MuonPFAnalyzer.h.

double MuonPFAnalyzer::theIsoCut [private]

Definition at line 101 of file MuonPFAnalyzer.h.

std::vector<std::string> MuonPFAnalyzer::theMuonKinds [private]

Definition at line 92 of file MuonPFAnalyzer.h.

std::map<std::string,std::map<std::string,MonitorElement*> > MuonPFAnalyzer::thePlots [private]

Definition at line 96 of file MuonPFAnalyzer.h.

Definition at line 97 of file MuonPFAnalyzer.h.

double MuonPFAnalyzer::theRecoGenR [private]

Definition at line 100 of file MuonPFAnalyzer.h.

Definition at line 88 of file MuonPFAnalyzer.h.

Definition at line 103 of file MuonPFAnalyzer.h.

Definition at line 89 of file MuonPFAnalyzer.h.