CMS 3D CMS Logo

GlobalMuonMatchAnalyzer.cc
Go to the documentation of this file.
1 
12 
13 // user include files
18 
21 
23 
27 
30 
31 #include <TH2.h>
32 
34 
35 {
36  iConfig = ps;
37  //now do what ever initialization is needed
40 
41  tkAssociatorToken_ = consumes<reco::TrackToTrackingParticleAssociator>(tkAssociatorName_);
42  muAssociatorToken_ = consumes<reco::TrackToTrackingParticleAssociator>(muAssociatorName_);
43 
44  subsystemname_ = iConfig.getUntrackedParameter<std::string>("subSystemFolder", "YourSubsystem");
49 
52 
53  tpToken_ = consumes<edm::View<reco::Track> >(tpName_);
54  tkToken_ = consumes<edm::View<reco::Track> >(tkName_);
55  staToken_ = consumes<edm::View<reco::Track> >(staName_);
56  glbToken_ = consumes<edm::View<reco::Track> >(glbName_);
57 }
58 
60  // do anything here that needs to be done at desctruction time
61  // (e.g. close files, deallocate resources etc.)
62 }
63 
64 //
65 // member functions
66 //
67 
68 // ------------ method called to for each event ------------
70  using namespace edm;
71  using namespace reco;
72 
74  iEvent.getByToken(tpToken_, tpHandle);
75  const TrackingParticleCollection tpColl = *(tpHandle.product());
76 
78  iEvent.getByToken(glbToken_, muHandle);
79  const reco::MuonTrackLinksCollection muColl = *(muHandle.product());
80 
81  Handle<View<Track> > staHandle;
82  iEvent.getByToken(staToken_, staHandle);
83 
84  Handle<View<Track> > glbHandle;
85  iEvent.getByToken(glbToken_, glbHandle);
86 
87  Handle<View<Track> > tkHandle;
88  iEvent.getByToken(tkToken_, tkHandle);
89 
91  iEvent.getByToken(tkAssociatorToken_, tkAssociator);
92 
93  // Mu Associator
95  iEvent.getByToken(muAssociatorToken_, muAssociator);
96 
97  reco::RecoToSimCollection tkrecoToSimCollection = tkAssociator->associateRecoToSim(tkHandle, tpHandle);
98  reco::SimToRecoCollection tksimToRecoCollection = tkAssociator->associateSimToReco(tkHandle, tpHandle);
99 
100  reco::RecoToSimCollection starecoToSimCollection = muAssociator->associateRecoToSim(staHandle, tpHandle);
101  reco::SimToRecoCollection stasimToRecoCollection = muAssociator->associateSimToReco(staHandle, tpHandle);
102 
103  reco::RecoToSimCollection glbrecoToSimCollection = muAssociator->associateRecoToSim(glbHandle, tpHandle);
104  reco::SimToRecoCollection glbsimToRecoCollection = muAssociator->associateSimToReco(glbHandle, tpHandle);
105 
106  for (TrackingParticleCollection::size_type i = 0; i < tpColl.size(); ++i) {
107  TrackingParticleRef tp(tpHandle, i);
108 
109  std::vector<std::pair<RefToBase<Track>, double> > rvGlb;
110  RefToBase<Track> rGlb;
111  if (glbsimToRecoCollection.find(tp) != glbsimToRecoCollection.end()) {
112  rvGlb = glbsimToRecoCollection[tp];
113  if (!rvGlb.empty()) {
114  rGlb = rvGlb.begin()->first;
115  }
116  }
117 
118  std::vector<std::pair<RefToBase<Track>, double> > rvSta;
119  RefToBase<Track> rSta;
120  if (stasimToRecoCollection.find(tp) != stasimToRecoCollection.end()) {
121  rvSta = stasimToRecoCollection[tp];
122  if (!rvSta.empty()) {
123  rSta = rvSta.begin()->first;
124  }
125  }
126 
127  std::vector<std::pair<RefToBase<Track>, double> > rvTk;
128  RefToBase<Track> rTk;
129  if (tksimToRecoCollection.find(tp) != tksimToRecoCollection.end()) {
130  rvTk = tksimToRecoCollection[tp];
131  if (!rvTk.empty()) {
132  rTk = rvTk.begin()->first;
133  }
134  }
135 
136  if (!rvSta.empty() && !rvTk.empty()) {
137  //should have matched
138  h_shouldMatch->Fill(rTk->eta(), rTk->pt());
139  }
140 
141  for (reco::MuonTrackLinksCollection::const_iterator links = muHandle->begin(); links != muHandle->end(); ++links) {
142  if (rGlb == RefToBase<Track>(links->globalTrack())) {
143  if (RefToBase<Track>(links->trackerTrack()) == rTk && RefToBase<Track>(links->standAloneTrack()) == rSta) {
144  //goodMatch
145  h_goodMatchSim->Fill(rGlb->eta(), rGlb->pt());
146  }
147  if (RefToBase<Track>(links->trackerTrack()) == rTk && RefToBase<Track>(links->standAloneTrack()) != rSta) {
148  //tkOnlyMatch
149  h_tkOnlySim->Fill(rGlb->eta(), rGlb->pt());
150  }
151  if (RefToBase<Track>(links->standAloneTrack()) == rSta && RefToBase<Track>(links->trackerTrack()) != rTk) {
152  //staOnlyMatch
153  h_staOnlySim->Fill(rGlb->eta(), rGlb->pt());
154  }
155  }
156  }
157  }
158 
160 
161  for (reco::MuonTrackLinksCollection::const_iterator links = muHandle->begin(); links != muHandle->end(); ++links) {
162  RefToBase<Track> glbRef = RefToBase<Track>(links->globalTrack());
163  RefToBase<Track> staRef = RefToBase<Track>(links->standAloneTrack());
164  RefToBase<Track> tkRef = RefToBase<Track>(links->trackerTrack());
165 
166  std::vector<std::pair<TrackingParticleRef, double> > tp1;
167  TrackingParticleRef tp1r;
168  if (glbrecoToSimCollection.find(glbRef) != glbrecoToSimCollection.end()) {
169  tp1 = glbrecoToSimCollection[glbRef];
170  if (!tp1.empty()) {
171  tp1r = tp1.begin()->first;
172  }
173  }
174 
175  std::vector<std::pair<TrackingParticleRef, double> > tp2;
176  TrackingParticleRef tp2r;
177  if (starecoToSimCollection.find(staRef) != starecoToSimCollection.end()) {
178  tp2 = starecoToSimCollection[staRef];
179  if (!tp2.empty()) {
180  tp2r = tp2.begin()->first;
181  }
182  }
183 
184  std::vector<std::pair<TrackingParticleRef, double> > tp3;
185  TrackingParticleRef tp3r;
186  if (tkrecoToSimCollection.find(tkRef) != tkrecoToSimCollection.end()) {
187  tp3 = tkrecoToSimCollection[tkRef];
188  if (!tp3.empty()) {
189  tp3r = tp3.begin()->first;
190  }
191  }
192 
193  if (!tp1.empty()) {
194  //was reconstructed
195  h_totReco->Fill(glbRef->eta(), glbRef->pt());
196  if (tp2r == tp3r) { // && tp1r == tp3r) {
197  //came from same TP
198  h_goodMatch->Fill(glbRef->eta(), glbRef->pt());
199  } else {
200  //mis-match
201  h_fakeMatch->Fill(glbRef->eta(), glbRef->pt());
202  }
203  }
204  }
205 }
206 
207 // ------------ method called once each job just after ending the event loop ------------
211 
214 
215  if (!out.empty() && dbe_)
216  dbe_->save(out);
217 }
218 
219 //void GlobalMuonMatchAnalyzer::beginRun(const edm::Run&, const edm::EventSetup& setup)
220 
222  edm::Run const &iRun,
223  edm::EventSetup const &iSetup) {
224  // Tk Associator
225 
226  ibooker.cd();
227  // Run histos only for dqmEndRun handling
228  ibooker.setScope(MonitorElementData::Scope::RUN);
229 
230  std::string dirName = "Matcher/";
231  // ibooker.setCurrentFolder("RecoMuonV/Matcher");
232  ibooker.setCurrentFolder(dirName);
233 
234  h_shouldMatch = ibooker.book2D("h_shouldMatch", "SIM associated to Tk and Sta", 50, -2.5, 2.5, 100, 0., 500.);
235  h_goodMatchSim = ibooker.book2D("h_goodMatchSim", "SIM associated to Glb Sta Tk", 50, -2.5, 2.5, 100, 0., 500.);
236  h_tkOnlySim = ibooker.book2D("h_tkOnlySim", "SIM associated to Glb Tk", 50, -2.5, 2.5, 100, 0., 500.);
237  h_staOnlySim = ibooker.book2D("h_staOnlySim", "SIM associated to Glb Sta", 50, -2.5, 2.5, 100, 0., 500.);
238 
239  h_totReco = ibooker.book2D("h_totReco", "Total Glb Reconstructed", 50, -2.5, 2.5, 100, 0., 500.);
240  h_goodMatch = ibooker.book2D("h_goodMatch", "Sta and Tk from same SIM", 50, -2.5, 2.5, 100, 0., 500.);
241  h_fakeMatch = ibooker.book2D("h_fakeMatch", "Sta and Tk not from same SIM", 50, -2.5, 2.5, 100, 0., 500.);
242 
243  h_effic = ibooker.book1D("h_effic", "Efficiency vs #eta", 50, -2.5, 2.5);
244  h_efficPt = ibooker.book1D("h_efficPt", "Efficiency vs p_{T}", 100, 0., 100.);
245 
246  h_fake = ibooker.book1D("h_fake", "Fake fraction vs #eta", 50, -2.5, 2.5);
247  h_fakePt = ibooker.book1D("h_fakePt", "Fake fraction vs p_{T}", 100, 0., 100.);
248 }
249 
251  MonitorElement *recoTH2,
252  MonitorElement *simTH2) {
253  TH2F *h1 = recoTH2->getTH2F();
254  TH1D *reco = h1->ProjectionX();
255 
256  TH2F *h2 = simTH2->getTH2F();
257  TH1D *sim = h2->ProjectionX();
258 
259  TH1F *hEff = (TH1F *)reco->Clone();
260 
261  hEff->Divide(sim);
262 
263  hEff->SetName("tmp_" + TString(reco->GetName()));
264 
265  // Set the error accordingly to binomial statistics
266  int nBinsEta = hEff->GetNbinsX();
267  for (int bin = 1; bin <= nBinsEta; bin++) {
268  float nSimHit = sim->GetBinContent(bin);
269  float eff = hEff->GetBinContent(bin);
270  float error = 0;
271  if (nSimHit != 0 && eff <= 1) {
272  error = sqrt(eff * (1 - eff) / nSimHit);
273  }
274  hEff->SetBinError(bin, error);
275  effHist->setBinContent(bin, eff);
276  effHist->setBinError(bin, error);
277  }
278 }
279 
281  MonitorElement *recoTH2,
282  MonitorElement *simTH2) {
283  TH2F *h1 = recoTH2->getTH2F();
284  TH1D *reco = h1->ProjectionY();
285 
286  TH2F *h2 = simTH2->getTH2F();
287  TH1D *sim = h2->ProjectionY();
288 
289  TH1F *hEff = (TH1F *)reco->Clone();
290 
291  hEff->Divide(sim);
292 
293  hEff->SetName("tmp_" + TString(reco->GetName()));
294 
295  // Set the error accordingly to binomial statistics
296  int nBinsPt = hEff->GetNbinsX();
297  for (int bin = 1; bin <= nBinsPt; bin++) {
298  float nSimHit = sim->GetBinContent(bin);
299  float eff = hEff->GetBinContent(bin);
300  float error = 0;
301  if (nSimHit != 0 && eff <= 1) {
302  error = sqrt(eff * (1 - eff) / nSimHit);
303  }
304  hEff->SetBinError(bin, error);
305  effHist->setBinContent(bin, eff);
306  effHist->setBinError(bin, error);
307  }
308 }
reco::SimToRecoCollection associateSimToReco(const edm::Handle< edm::View< reco::Track >> &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > muAssociatorToken_
reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track >> &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
T const * product() const
Definition: Handle.h:70
virtual MonitorElementData::Scope setScope(MonitorElementData::Scope newscope)
Definition: DQMStore.cc:50
uint16_t size_type
void dqmEndRun(const edm::Run &, const edm::EventSetup &) override
const_iterator find(const key_type &k) const
find element with specified reference key
std::vector< MuonTrackLinks > MuonTrackLinksCollection
collection of MuonTrackLinks
Definition: MuonFwd.h:22
const_iterator end() const
last iterator over the map (read only)
T getUntrackedParameter(std::string const &, T const &) const
void Fill(long long x)
double pt() const
track transverse momentum
Definition: TrackBase.h:637
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< edm::View< reco::Track > > glbToken_
T sqrt(T t)
Definition: SSEVec.h:23
void computeEfficiencyPt(MonitorElement *, MonitorElement *recoTH2, MonitorElement *simTH2)
edm::EDGetTokenT< edm::View< reco::Track > > staToken_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
void analyze(const edm::Event &, const edm::EventSetup &) override
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
GlobalMuonMatchAnalyzer(const edm::ParameterSet &)
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
DQM_DEPRECATED void save(std::string const &filename, std::string const &path="")
Definition: DQMStore.cc:824
fixed size matrix
HLT enums.
edm::EDGetTokenT< edm::View< reco::Track > > tpToken_
virtual void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
std::vector< TrackingParticle > TrackingParticleCollection
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
void computeEfficiencyEta(MonitorElement *, MonitorElement *recoTH2, MonitorElement *simTH2)
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > tkAssociatorToken_
Definition: Run.h:45
edm::EDGetTokenT< edm::View< reco::Track > > tkToken_