CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GlobalMuonMatchAnalyzer.cc
Go to the documentation of this file.
1 
12 
13 // user include files
19 
22 
26 
30 
31 
34 
35 
36 #include <TH2.h>
37 
38 
40 
41 {
42  iConfig = ps;
43  //now do what ever initialization is needed
46 
47 subsystemname_ = iConfig.getUntrackedParameter<std::string>("subSystemFolder", "YourSubsystem") ;
52 
55 
56  tpToken_ = consumes<edm::View<reco::Track> >(tpName_);
57  tkToken_ = consumes<edm::View<reco::Track> >(tkName_);
58  staToken_ = consumes<edm::View<reco::Track> >(staName_);
59  glbToken_ = consumes<edm::View<reco::Track> >(glbName_);
60 }
61 
62 
64 {
65 
66  // do anything here that needs to be done at desctruction time
67  // (e.g. close files, deallocate resources etc.)
68 
69 }
70 
71 
72 //
73 // member functions
74 //
75 
76 // ------------ method called to for each event ------------
77 void
79 {
80  using namespace edm;
81  using namespace reco;
82 
84  iEvent.getByToken(tpToken_,tpHandle);
85  const TrackingParticleCollection tpColl = *(tpHandle.product());
86 
88  iEvent.getByToken(glbToken_,muHandle);
89  const reco::MuonTrackLinksCollection muColl = *(muHandle.product());
90 
91  Handle<View<Track> > staHandle;
92  iEvent.getByToken(staToken_,staHandle);
93 
94  Handle<View<Track> > glbHandle;
95  iEvent.getByToken(glbToken_,glbHandle);
96 
97  Handle<View<Track> > tkHandle;
98  iEvent.getByToken(tkToken_,tkHandle);
99 
100  reco::RecoToSimCollection tkrecoToSimCollection;
101  reco::SimToRecoCollection tksimToRecoCollection;
102  tkrecoToSimCollection = tkAssociator_->associateRecoToSim(tkHandle,tpHandle,&iEvent,&iSetup);
103  tksimToRecoCollection = tkAssociator_->associateSimToReco(tkHandle,tpHandle,&iEvent,&iSetup);
104 
105  reco::RecoToSimCollection starecoToSimCollection;
106  reco::SimToRecoCollection stasimToRecoCollection;
107  starecoToSimCollection = muAssociator_->associateRecoToSim(staHandle,tpHandle,&iEvent,&iSetup);
108  stasimToRecoCollection = muAssociator_->associateSimToReco(staHandle,tpHandle,&iEvent,&iSetup);
109 
110  reco::RecoToSimCollection glbrecoToSimCollection;
111  reco::SimToRecoCollection glbsimToRecoCollection;
112  glbrecoToSimCollection = muAssociator_->associateRecoToSim(glbHandle,tpHandle,&iEvent,&iSetup);
113  glbsimToRecoCollection = muAssociator_->associateSimToReco(glbHandle,tpHandle,&iEvent,&iSetup);
114 
115 
116  for (TrackingParticleCollection::size_type i=0; i<tpColl.size(); ++i){
117  TrackingParticleRef tp(tpHandle,i);
118 
119  std::vector<std::pair<RefToBase<Track>, double> > rvGlb;
120  RefToBase<Track> rGlb;
121  if(glbsimToRecoCollection.find(tp) != glbsimToRecoCollection.end()){
122  rvGlb = glbsimToRecoCollection[tp];
123  if(rvGlb.size() != 0) {
124  rGlb = rvGlb.begin()->first;
125  }
126  }
127 
128  std::vector<std::pair<RefToBase<Track>, double> > rvSta;
129  RefToBase<Track> rSta;
130  if(stasimToRecoCollection.find(tp) != stasimToRecoCollection.end()){
131  rvSta = stasimToRecoCollection[tp];
132  if(rvSta.size() != 0) {
133  rSta = rvSta.begin()->first;
134  }
135  }
136 
137  std::vector<std::pair<RefToBase<Track>, double> > rvTk;
138  RefToBase<Track> rTk;
139  if(tksimToRecoCollection.find(tp) != tksimToRecoCollection.end()){
140  rvTk = tksimToRecoCollection[tp];
141  if(rvTk.size() != 0) {
142  rTk = rvTk.begin()->first;
143  }
144  }
145 
146  if( rvSta.size() != 0 && rvTk.size() != 0 ){
147  //should have matched
148  h_shouldMatch->Fill(rTk->eta(),rTk->pt());
149  }
150 
151  for ( reco::MuonTrackLinksCollection::const_iterator links = muHandle->begin(); links != muHandle->end(); ++links ) {
152  if( rGlb == RefToBase<Track>(links->globalTrack() ) ) {
153  if( RefToBase<Track>(links->trackerTrack() ) == rTk &&
154  RefToBase<Track>(links->standAloneTrack() ) == rSta ) {
155  //goodMatch
156  h_goodMatchSim->Fill(rGlb->eta(),rGlb->pt());
157  }
158  if ( RefToBase<Track>(links->trackerTrack() ) == rTk &&
159  RefToBase<Track>(links->standAloneTrack() ) != rSta ) {
160  //tkOnlyMatch
161  h_tkOnlySim->Fill(rGlb->eta(),rGlb->pt());
162  }
163  if ( RefToBase<Track>(links->standAloneTrack() ) == rSta &&
164  RefToBase<Track>(links->trackerTrack() ) != rTk ) {
165  //staOnlyMatch
166  h_staOnlySim->Fill(rGlb->eta(),rGlb->pt());
167  }
168  }
169  }
170 
171  }
172 
174 
175  for ( reco::MuonTrackLinksCollection::const_iterator links = muHandle->begin(); links != muHandle->end(); ++links ) {
176  RefToBase<Track> glbRef = RefToBase<Track>(links->globalTrack() );
177  RefToBase<Track> staRef = RefToBase<Track>(links->standAloneTrack() );
178  RefToBase<Track> tkRef = RefToBase<Track>(links->trackerTrack() );
179 
180  std::vector<std::pair<TrackingParticleRef, double> > tp1;
181  TrackingParticleRef tp1r;
182  if(glbrecoToSimCollection.find(glbRef) != glbrecoToSimCollection.end()){
183  tp1 = glbrecoToSimCollection[glbRef];
184  if(tp1.size() != 0) {
185  tp1r = tp1.begin()->first;
186  }
187  }
188 
189  std::vector<std::pair<TrackingParticleRef, double> > tp2;
190  TrackingParticleRef tp2r;
191  if(starecoToSimCollection.find(staRef) != starecoToSimCollection.end()){
192  tp2 = starecoToSimCollection[staRef];
193  if(tp2.size() != 0) {
194  tp2r = tp2.begin()->first;
195  }
196  }
197 
198  std::vector<std::pair<TrackingParticleRef, double> > tp3;
199  TrackingParticleRef tp3r;
200  if(tkrecoToSimCollection.find(tkRef) != tkrecoToSimCollection.end()){
201  tp3 = tkrecoToSimCollection[tkRef];
202  if(tp3.size() != 0) {
203  tp3r = tp3.begin()->first;
204  }
205  }
206 
207 
208  if(tp1.size() != 0) {
209  //was reconstructed
210  h_totReco->Fill(glbRef->eta(),glbRef->pt());
211  if(tp2r == tp3r) { // && tp1r == tp3r) {
212  //came from same TP
213  h_goodMatch->Fill(glbRef->eta(),glbRef->pt());
214  } else {
215  //mis-match
216  h_fakeMatch->Fill(glbRef->eta(),glbRef->pt());
217  }
218  }
219 
220  }
221 
222 
223 }
224 
225 
226 // ------------ method called once each job just before starting event loop ------------
227 void
229 {
230 }
232 // ------------ method called once each job just after ending the event loop ------------
236 
239 
240  if( out.size() != 0 && dbe_ ) dbe_->save(out);
241 }
242 
243 //void GlobalMuonMatchAnalyzer::beginRun(const edm::Run&, const edm::EventSetup& setup)
244 
246  edm::Run const & iRun,
247  edm::EventSetup const & iSetup )
248 {
249  // Tk Associator
250 
251  edm::ESHandle<TrackAssociatorBase> tkassociatorHandle;
252  iSetup.get<TrackAssociatorRecord>().get(tkAssociatorName_,tkassociatorHandle);
253  tkAssociator_ = tkassociatorHandle.product();
254 
255  // Mu Associator
256  edm::ESHandle<TrackAssociatorBase> muassociatorHandle;
257  iSetup.get<TrackAssociatorRecord>().get(muAssociatorName_,muassociatorHandle);
258  muAssociator_ = muassociatorHandle.product();
259  ibooker.cd();
260  std::string dirName="Matcher/";
261  // ibooker.setCurrentFolder("RecoMuonV/Matcher");
262  ibooker.setCurrentFolder(dirName.c_str()) ;
263 
264  h_shouldMatch = ibooker.book2D("h_shouldMatch","SIM associated to Tk and Sta",50,-2.5,2.5,100,0.,500.);
265  h_goodMatchSim = ibooker.book2D("h_goodMatchSim","SIM associated to Glb Sta Tk",50,-2.5,2.5,100,0.,500.);
266  h_tkOnlySim = ibooker.book2D("h_tkOnlySim","SIM associated to Glb Tk",50,-2.5,2.5,100,0.,500.);
267  h_staOnlySim = ibooker.book2D("h_staOnlySim","SIM associated to Glb Sta",50,-2.5,2.5,100,0.,500.);
268 
269  h_totReco = ibooker.book2D("h_totReco","Total Glb Reconstructed",50,-2.5,2.5,100,0.,500.);
270  h_goodMatch = ibooker.book2D("h_goodMatch","Sta and Tk from same SIM",50,-2.5,2.5,100, 0., 500.);
271  h_fakeMatch = ibooker.book2D("h_fakeMatch","Sta and Tk not from same SIM",50,-2.5,2.5,100,0.,500.);
272 
273  h_effic = ibooker.book1D("h_effic","Efficiency vs #eta",50,-2.5,2.5);
274  h_efficPt = ibooker.book1D("h_efficPt","Efficiency vs p_{T}",100,0.,100.);
275 
276  h_fake = ibooker.book1D("h_fake","Fake fraction vs #eta",50,-2.5,2.5);
277  h_fakePt = ibooker.book1D("h_fakePt","Fake fraction vs p_{T}",100,0.,100.);
278 
279 }
280 
281 
282 
284  TH2F * h1 = recoTH2->getTH2F();
285  TH1D* reco = h1->ProjectionX();
286 
287  TH2F * h2 =simTH2->getTH2F();
288  TH1D* sim = h2 ->ProjectionX();
289 
290 
291  TH1F *hEff = (TH1F*) reco->Clone();
292 
293  hEff->Divide(sim);
294 
295  hEff->SetName("tmp_"+TString(reco->GetName()));
296 
297  // Set the error accordingly to binomial statistics
298  int nBinsEta = hEff->GetNbinsX();
299  for(int bin = 1; bin <= nBinsEta; bin++) {
300  float nSimHit = sim->GetBinContent(bin);
301  float eff = hEff->GetBinContent(bin);
302  float error = 0;
303  if(nSimHit != 0 && eff <= 1) {
304  error = sqrt(eff*(1-eff)/nSimHit);
305  }
306  hEff->SetBinError(bin, error);
307  effHist->setBinContent(bin,eff);
308  effHist->setBinError(bin,error);
309  }
310 
311 }
312 
314  TH2F * h1 = recoTH2->getTH2F();
315  TH1D* reco = h1->ProjectionY();
316 
317  TH2F * h2 = simTH2->getTH2F();
318  TH1D* sim = h2 ->ProjectionY();
319 
320 
321  TH1F *hEff = (TH1F*) reco->Clone();
322 
323  hEff->Divide(sim);
324 
325  hEff->SetName("tmp_"+TString(reco->GetName()));
326 
327  // Set the error accordingly to binomial statistics
328  int nBinsPt = hEff->GetNbinsX();
329  for(int bin = 1; bin <= nBinsPt; bin++) {
330  float nSimHit = sim->GetBinContent(bin);
331  float eff = hEff->GetBinContent(bin);
332  float error = 0;
333  if(nSimHit != 0 && eff <= 1) {
334  error = sqrt(eff*(1-eff)/nSimHit);
335  }
336  hEff->SetBinError(bin, error);
337  effHist->setBinContent(bin,eff);
338  effHist->setBinError(bin,error);
339  }
340 
341 }
342 
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
const TrackAssociatorBase * muAssociator_
void setBinContent(int binx, double content)
set content of bin (1-D)
std::vector< TrackingParticle > TrackingParticleCollection
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void cd(void)
Definition: DQMStore.cc:266
Definition: sim.h:19
uint16_t size_type
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::vector< MuonTrackLinks > MuonTrackLinksCollection
collection of MuonTrackLinks
Definition: MuonFwd.h:22
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:139
edm::EDGetTokenT< edm::View< reco::Track > > glbToken_
T sqrt(T t)
Definition: SSEVec.h:48
void computeEfficiencyPt(MonitorElement *, MonitorElement *recoTH2, MonitorElement *simTH2)
double pt() const
track transverse momentum
Definition: TrackBase.h:129
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:113
void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
edm::EDGetTokenT< edm::View< reco::Track > > staToken_
virtual reco::RecoToSimCollection associateRecoToSim(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event, const edm::EventSetup *setup) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
virtual reco::SimToRecoCollection associateSimToReco(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event, const edm::EventSetup *setup) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:131
GlobalMuonMatchAnalyzer(const edm::ParameterSet &)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
const TrackAssociatorBase * tkAssociator_
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, const uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE", const bool resetMEsAfterWriting=false)
Definition: DQMStore.cc:2490
edm::EDGetTokenT< edm::View< reco::Track > > tpToken_
TH2F * getTH2F(void) const
void computeEfficiencyEta(MonitorElement *, MonitorElement *recoTH2, MonitorElement *simTH2)
Definition: Run.h:41
edm::EDGetTokenT< edm::View< reco::Track > > tkToken_