CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTTauDQMTrkPlotter.cc
Go to the documentation of this file.
2 
3 HLTTauDQMTrkPlotter::HLTTauDQMTrkPlotter(const edm::ParameterSet& iConfig, int etbins, int etabins, int phibins, double maxpt, bool ref, double dr, std::string dqmBaseFolder ) {
4  //Initialize Plotter
5  name_ = "HLTTauDQMTrkPlotter";
6 
7  //Process PSet
8  try {
9  jetTagSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("ConeIsolation");
10  isolJets_ = iConfig.getUntrackedParameter<edm::InputTag>("IsolatedJets");
11  triggerTag_ = iConfig.getUntrackedParameter<std::string>("DQMFolder");
12  triggerTagAlias_ = iConfig.getUntrackedParameter<std::string>("Alias","");
13  type_ = iConfig.getUntrackedParameter<std::string>("Type");
14  mcMatch_ = dr;
15  EtMax_ = maxpt;
16  NPtBins_ = etbins;
17  NEtaBins_ = etabins;
18  NPhiBins_ = phibins;
19  dqmBaseFolder_ = dqmBaseFolder;
20  doRef_ = ref;
21  validity_ = true;
22  } catch ( cms::Exception &e ) {
23  edm::LogWarning("HLTTauDQMTrkPlotter::HLTTauDQMTrkPlotter") << e.what() << std::endl;
24  validity_ = false;
25  return;
26  }
27 
28  if (store_) {
29  //Create the histograms
32 
33  jetEt = store_->book1D((type_+"TauEt").c_str(), "#tau E_{T}",NPtBins_,0,EtMax_);
34  jetEta = store_->book1D((type_+"TauEta").c_str(), "#tau #eta", NEtaBins_, -2.5, 2.5);
35  jetPhi = store_->book1D((type_+"TauPhi").c_str(), "#tau #phi", NPhiBins_, -3.2, 3.2);
36  isoJetEt = store_->book1D((type_+"IsolJetEt").c_str(), "Selected Jet E_{T}", NPtBins_, 0,EtMax_);
37  isoJetEta = store_->book1D((type_+"IsolJetEta").c_str(), "Selected Jet #eta", NEtaBins_, -2.5, 2.5);
38  isoJetPhi = store_->book1D((type_+"IsolJetPhi").c_str(), "Selected jet #phi", NPhiBins_, -3.2, 3.2);
39 
40  nPxlTrksInL25Jet = store_->book1D((type_+"nTracks").c_str(), "# RECO Tracks", 30, 0, 30);
41  nQPxlTrksInL25Jet = store_->book1D((type_+"nQTracks").c_str(),"# Quality RECO Tracks", 15, 0, 15);
42  signalLeadTrkPt = store_->book1D((type_+"LeadTrackPt").c_str(), "Lead Track p_{T}", 75, 0, 150);
43  hasLeadTrack = store_->book1D((type_+"HasLeadTrack").c_str(), "Lead Track ?", 2, 0, 2);
44 
45  EtEffNum = store_->book1D((type_+"TauEtEffNum").c_str(),"Efficiency vs E_{T} (Numerator)",NPtBins_,0,EtMax_);
46  EtEffNum->getTH1F()->Sumw2();
47 
48  EtEffDenom = store_->book1D((type_+"TauEtEffDenom").c_str(),"Efficiency vs E_{T} (Denominator)",NPtBins_,0,EtMax_);
49  EtEffDenom->getTH1F()->Sumw2();
50 
51  EtaEffNum = store_->book1D((type_+"TauEtaEffNum").c_str(),"Efficiency vs #eta (Numerator)",NEtaBins_,-2.5,2.5);
52  EtaEffNum->getTH1F()->Sumw2();
53 
54  EtaEffDenom = store_->book1D((type_+"TauEtaEffDenom").c_str(),"Efficiency vs #eta (Denominator)",NEtaBins_,-2.5,2.5);
55  EtaEffDenom->getTH1F()->Sumw2();
56 
57  PhiEffNum = store_->book1D((type_+"TauPhiEffNum").c_str(),"Efficiency vs #phi (Numerator)",NPhiBins_,-3.2,3.2);
58  PhiEffNum->getTH1F()->Sumw2();
59 
60  PhiEffDenom = store_->book1D((type_+"TauPhiEffDenom").c_str(),"Efficiency vs #phi (Denominator)",NPhiBins_,-3.2,3.2);
61  PhiEffDenom->getTH1F()->Sumw2();
62  }
63 }
64 
65 
67 }
68 
69 void HLTTauDQMTrkPlotter::analyze( const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::map<int,LVColl>& mcInfo ) {
70  using namespace edm;
71  using namespace reco;
72 
73  //Tau reference
74  std::map<int,LVColl>::const_iterator iref;
75  iref = mcInfo.find(15);
76 
78  Handle<CaloJetCollection> isolJets;
79 
80  bool gotL3 = iEvent.getByLabel(jetTagSrc_, tauTagInfos) && tauTagInfos.isValid();
81 
82  if ( gotL3 ) {
83  for ( unsigned int i = 0; i < tauTagInfos->size(); ++i ) {
84  IsolatedTauTagInfo tauTagInfo = (*tauTagInfos)[i];
85  if ( &(*tauTagInfo.jet()) ) {
86  LV theJet = tauTagInfo.jet()->p4();
87 
88  std::pair<bool,LV> m(false,LV());
89  if ( iref != mcInfo.end() ) m = match(theJet,iref->second,mcMatch_);
90 
91  if ( (doRef_ && m.first) || (!doRef_) ) {
92  jetEt->Fill(theJet.Et());
93  jetEta->Fill(theJet.Eta());
94  jetPhi->Fill(theJet.Phi());
95  nPxlTrksInL25Jet->Fill(tauTagInfo.allTracks().size());
96  nQPxlTrksInL25Jet->Fill(tauTagInfo.selectedTracks().size());
97 
98  const TrackRef leadTrk = tauTagInfo.leadingSignalTrack();
99  if ( !leadTrk ) {
100  hasLeadTrack->Fill(0.5);
101  } else {
102  hasLeadTrack->Fill(1.5);
103  signalLeadTrkPt->Fill(leadTrk->pt());
104  }
105 
106  LV refV;
107  if ( doRef_ ) {
108  refV = m.second;
109  } else {
110  refV = theJet;
111  }
112 
113  EtEffDenom->Fill(refV.pt());
114  EtaEffDenom->Fill(refV.eta());
115  PhiEffDenom->Fill(refV.phi());
116 
117  bool gotIsoL3 = iEvent.getByLabel(isolJets_, isolJets) && isolJets.isValid();
118 
119  if ( gotIsoL3 ) {
120  if ( matchJet(*(tauTagInfo.jet()),*isolJets) ) {
121  isoJetEta->Fill(theJet.Eta());
122  isoJetEt->Fill(theJet.Et());
123  isoJetPhi->Fill(theJet.Phi());
124 
125  EtEffNum->Fill(refV.pt());
126  EtaEffNum->Fill(refV.eta());
127  PhiEffNum->Fill(refV.phi());
128  }
129  }
130  }
131  }
132  }
133  }
134 }
135 
137  //Loop On the Collection and see if your tau jet is matched to one there
138  //Also find the nearest Matched MC Particle to your Jet (to be complete)
139 
140  bool matched = false;
141 
142  for ( reco::CaloJetCollection::const_iterator it = McInfo.begin(); it != McInfo.end(); ++it ) {
143  double delta = ROOT::Math::VectorUtil::DeltaR(jet.p4().Vect(),it->p4().Vect());
144  if ( delta < mcMatch_ ) {
145  matched = true;
146  }
147  }
148  return matched;
149 }
virtual char const * what() const
Definition: Exception.cc:141
dbl * delta
Definition: mlp_gen.cc:36
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
MonitorElement * EtaEffNum
HLTTauDQMTrkPlotter(const edm::ParameterSet &, int, int, int, double, bool, double, std::string)
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:722
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
Base class for all types of Jets.
Definition: Jet.h:21
MonitorElement * jetPhi
MonitorElement * jetEt
void Fill(long long x)
const TrackRefVector selectedTracks() const
MonitorElement * PhiEffDenom
MonitorElement * isoJetEta
math::XYZTLorentzVectorD LV
int iEvent
Definition: GenABIO.cc:243
void removeContents(void)
erase all monitoring elements in current directory (not including subfolders);
Definition: DQMStore.cc:2569
bool matchJet(const reco::Jet &, const reco::CaloJetCollection &)
std::string triggerTagAlias_
MonitorElement * jetEta
MonitorElement * nPxlTrksInL25Jet
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
virtual edm::RefToBase< Jet > jet(void) const
returns a polymorphic reference to the tagged jet
Definition: JTATagInfo.h:20
MonitorElement * EtEffNum
MonitorElement * nQPxlTrksInL25Jet
MonitorElement * isoJetEt
MonitorElement * hasLeadTrack
MonitorElement * signalLeadTrkPt
std::string triggerTag()
void analyze(const edm::Event &, const edm::EventSetup &, const std::map< int, LVColl > &)
MonitorElement * PhiEffNum
TH1F * getTH1F(void) const
const TrackRef leadingSignalTrack() const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
MonitorElement * EtEffDenom
std::pair< bool, LV > match(const LV &, const LVColl &, double)
MonitorElement * EtaEffDenom
std::string triggerTag_
const TrackRefVector allTracks() const
MonitorElement * isoJetPhi
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:434
std::string dqmBaseFolder_
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects