CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

L25TauAnalyzer Class Reference

#include <HLTriggerOffline/Tau/src/L25TauAnalyzer.cc>

Inheritance diagram for L25TauAnalyzer:
edm::EDAnalyzer

List of all members.

Public Member Functions

 L25TauAnalyzer (const edm::ParameterSet &)
 ~L25TauAnalyzer ()

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void endJob ()
MatchElementL25 match (const reco::Jet &, const LVColl &)
float trackDrRMS (const reco::IsolatedTauTagInfo &, const reco::TrackRefVector &)

Private Attributes

float emf
bool hasLeadTrk
float isolationCone_
float jetEt
float jetEta
float jetMCEt
float jetMCEta
edm::InputTag jetMCTagSrc_
edm::InputTag jetTagSrc_
TFile * l25file
TTree * l25tree
float leadSignalTrackPt
float leadTrkJetDeltaR
float minTrackPt_
int numPixTrkInJet
int numQPixTrkInAnnulus
int numQPixTrkInJet
int numQPixTrkInSignalCone
std::string rootFile_
bool signal_
float signalCone_
float trkDrRMS
float trkDrRMSA

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 63 of file L25TauAnalyzer.h.


Constructor & Destructor Documentation

L25TauAnalyzer::L25TauAnalyzer ( const edm::ParameterSet iConfig) [explicit]

Definition at line 30 of file L25TauAnalyzer.cc.

References edm::ParameterSet::getParameter(), and reco::btau::jetEta.

                                                            {
  jetTagSrc_ = iConfig.getParameter<InputTag>("JetTagProd");
  jetMCTagSrc_ = iConfig.getParameter<InputTag>("JetMCTagProd");
  rootFile_ = iConfig.getParameter<string>("outputFileName");
  signal_ = iConfig.getParameter<bool>("Signal");
  minTrackPt_ =iConfig.getParameter<double>("MinTrackPt");
  signalCone_ =iConfig.getParameter<double>("SignalCone");
  isolationCone_ =iConfig.getParameter<double>("IsolationCone");
 
  l25file = new TFile(rootFile_.c_str(),"recreate");
  l25tree = new TTree("l25tree","Level 2.5 Tau Tree");

  jetEt=0.;
  jetEta=0.;
  jetMCEt=0.;
  jetMCEta=0.;
  leadSignalTrackPt=0.;
  trkDrRMS=0.;
  trkDrRMSA=0.;
  leadTrkJetDeltaR=0.;
  numPixTrkInJet=0;
  numQPixTrkInJet=0;
  numQPixTrkInSignalCone=0;
  numQPixTrkInAnnulus=0;
  hasLeadTrk=0;
  emf=0;


  l25tree->Branch("jetEt", &jetEt, "jetEt/F");
  l25tree->Branch("jetEMF", &emf, "jetEMF/F");
  l25tree->Branch("jetEta", &jetEta, "jetEta/F");
  l25tree->Branch("jetMCEt", &jetMCEt, "jetMCEt/F");
  l25tree->Branch("jetMCEta", &jetMCEta, "jetMCEta/F");
  l25tree->Branch("leadTrackPt", &leadSignalTrackPt, "leadTrackPt/F");
  l25tree->Branch("trackDeltaRRMS", &trkDrRMS, "trackDeltaRRMS/F");
  l25tree->Branch("trackDeltaRRMSAll", &trkDrRMSA, "trackDeltaRRMSAll/F");
  l25tree->Branch("matchingCone", &leadTrkJetDeltaR, "matchingCone/F");
  l25tree->Branch("nTracks", &numPixTrkInJet, "nTracks/I");
  l25tree->Branch("nQTracks", &numQPixTrkInJet, "nQTracks/I");
  l25tree->Branch("nQTracksInSignal", &numQPixTrkInSignalCone, "nQTracksInSignal/I");
  l25tree->Branch("nQTracksInAnnulus", &numQPixTrkInAnnulus, "nQTracksInAnnulus/I");
  l25tree->Branch("hasLeadTrack", &hasLeadTrk, "hasLeadTrack/B");
}
L25TauAnalyzer::~L25TauAnalyzer ( )

Definition at line 75 of file L25TauAnalyzer.cc.

                               {

}

Member Function Documentation

void L25TauAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 82 of file L25TauAnalyzer.cc.

References MatchElementL25::deltar, reco::CaloJet::emEnergyFraction(), edm::Event::getByLabel(), i, reco::btau::jetEta, m, edm::match(), MatchElementL25::matched, MatchElementL25::mcEt, MatchElementL25::mcEta, and dt_dqm_sourceclient_common_cff::reco.

                                                                               {


   Handle<IsolatedTauTagInfoCollection> tauTagInfo;
   using namespace reco;
   Handle<LVColl> mcInfo;                                                        

   if(signal_){                                                                  
      iEvent.getByLabel(jetMCTagSrc_, mcInfo);                           
   }                                                                     
   

   //get isolation tag
   if(iEvent.getByLabel(jetTagSrc_, tauTagInfo))
     {
       if(tauTagInfo->size()>0)
         for(IsolatedTauTagInfoCollection::const_iterator i = tauTagInfo->begin();i!=tauTagInfo->end();++i)
           {
             MatchElementL25 m;
             m.matched=false;
             m.mcEt=0;
             m.mcEta=0;
             m.deltar=0;

             if(signal_)
               {
                 if(iEvent.getByLabel(jetMCTagSrc_,mcInfo))
                   m=match(*(i->jet()),*mcInfo);
               }

             if((signal_&&m.matched)||(!signal_))
               {

                 const Jet* Jet =i->jet().get();
                 const CaloJet* calojet = dynamic_cast<const CaloJet*>(Jet);
                 emf = calojet->emEnergyFraction();

                 
                 jetEt=i->jet()->et();
                 jetEta=i->jet()->eta();
                 jetMCEt=m.mcEt;
                 jetMCEta=m.mcEta;
                 numPixTrkInJet = i->allTracks().size();
                 numQPixTrkInJet = i->selectedTracks().size();
                 trkDrRMS =trackDrRMS(*i,i->selectedTracks());
                 trkDrRMSA =trackDrRMS(*i,i->allTracks());
 

                 //Search Leading Track
                 const TrackRef leadTk = i->leadingSignalTrack();
                 if(!leadTk)
                   {
                     numQPixTrkInSignalCone=0;
                     numQPixTrkInAnnulus=0;
                     leadSignalTrackPt=0;
                     leadTrkJetDeltaR=0;
                     hasLeadTrk=false;
                   }
                 else
                   {
                     leadTrkJetDeltaR = ROOT::Math::VectorUtil::DeltaR(i->jet()->p4().Vect(), leadTk->momentum());                            
                     leadSignalTrackPt=leadTk->pt();
                     numQPixTrkInSignalCone=(i->tracksInCone(leadTk->momentum(), signalCone_, minTrackPt_)).size();
                     numQPixTrkInAnnulus=(i->tracksInCone(leadTk->momentum(), isolationCone_, minTrackPt_)).size()-numQPixTrkInSignalCone;
                     hasLeadTrk=true;
                   }
               }
              l25tree->Fill();
           }
     }
}
void L25TauAnalyzer::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 234 of file L25TauAnalyzer.cc.

                              {
}
void L25TauAnalyzer::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 238 of file L25TauAnalyzer.cc.

                            {
   l25file->Write();
}
MatchElementL25 L25TauAnalyzer::match ( const reco::Jet jet,
const LVColl McInfo 
) [private]

Definition at line 193 of file L25TauAnalyzer.cc.

References delta, MatchElementL25::deltar, edm::match(), MatchElementL25::matched, MatchElementL25::mcEt, MatchElementL25::mcEta, and reco::LeafCandidate::p4().

{

  //Loop On the Collection and see if your tau jet is matched to one there
 //Also find the nearest Matched MC Particle to your Jet (to be complete)
 
 bool matched=false;
 double delta_min=100.;
 double mceta=0;
 double mcet=0;
 
 double matchDr=0.3;

 if(McInfo.size()>0)
  for(std::vector<LV>::const_iterator it = McInfo.begin();it!=McInfo.end();++it)
   {
          double delta = ROOT::Math::VectorUtil::DeltaR(jet.p4().Vect(),*it);
          if(delta<matchDr)
            {
              matched=true;
              if(delta<delta_min)
                {
                  delta_min=delta;
                  mceta=it->eta();
                  mcet=it->Et();
                }
            }
   }

  //Create Struct and send it out!
  MatchElementL25 match;
  match.matched=matched;
  match.deltar=delta_min;
  match.mcEta = mceta;
  match.mcEt = mcet;

 return match;
}
float L25TauAnalyzer::trackDrRMS ( const reco::IsolatedTauTagInfo info,
const reco::TrackRefVector tracks 
) [private]

Definition at line 158 of file L25TauAnalyzer.cc.

References i, funct::pow(), plotscripts::rms(), and edm::RefVector< C, T, F >::size().

{
  float rms = 0.;
  float sumet = 0.;
  
  
  //First find the weighted track
  if(tracks.size()>0)
    {
      math::XYZVector center = tracks[0]->momentum();

      for(size_t i = 1;i<tracks.size();++i)
        {
          center+=tracks[i]->momentum();
        }

      //Now calculate DeltaR

      for(size_t i = 0;i<tracks.size();++i)
        {
          rms+= tracks[i]->pt()*pow(ROOT::Math::VectorUtil::DeltaR(tracks[i]->momentum(),center),2);
          sumet+=tracks[i]->pt();
        }
      
    }
     
  if(sumet==0.)
    sumet=1;
  
  return rms/sumet;
}

Member Data Documentation

float L25TauAnalyzer::emf [private]

Definition at line 100 of file L25TauAnalyzer.h.

Definition at line 101 of file L25TauAnalyzer.h.

Definition at line 82 of file L25TauAnalyzer.h.

float L25TauAnalyzer::jetEt [private]

Definition at line 92 of file L25TauAnalyzer.h.

float L25TauAnalyzer::jetEta [private]

Definition at line 93 of file L25TauAnalyzer.h.

float L25TauAnalyzer::jetMCEt [private]

Definition at line 94 of file L25TauAnalyzer.h.

float L25TauAnalyzer::jetMCEta [private]

Definition at line 95 of file L25TauAnalyzer.h.

Definition at line 77 of file L25TauAnalyzer.h.

Definition at line 76 of file L25TauAnalyzer.h.

TFile* L25TauAnalyzer::l25file [private]

Definition at line 85 of file L25TauAnalyzer.h.

TTree* L25TauAnalyzer::l25tree [private]

Definition at line 86 of file L25TauAnalyzer.h.

Definition at line 98 of file L25TauAnalyzer.h.

Definition at line 99 of file L25TauAnalyzer.h.

float L25TauAnalyzer::minTrackPt_ [private]

Definition at line 80 of file L25TauAnalyzer.h.

Definition at line 88 of file L25TauAnalyzer.h.

Definition at line 91 of file L25TauAnalyzer.h.

Definition at line 89 of file L25TauAnalyzer.h.

Definition at line 90 of file L25TauAnalyzer.h.

std::string L25TauAnalyzer::rootFile_ [private]

Definition at line 78 of file L25TauAnalyzer.h.

bool L25TauAnalyzer::signal_ [private]

Definition at line 79 of file L25TauAnalyzer.h.

float L25TauAnalyzer::signalCone_ [private]

Definition at line 81 of file L25TauAnalyzer.h.

float L25TauAnalyzer::trkDrRMS [private]

Definition at line 96 of file L25TauAnalyzer.h.

float L25TauAnalyzer::trkDrRMSA [private]

Definition at line 97 of file L25TauAnalyzer.h.