CMS 3D CMS Logo

OverlapProblemTPAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: OverlapProblemTPAnalyzer
4 // Class: OverlapProblemTPAnalyzer
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Thu Dec 16 16:32:56 CEST 2010
16 // $Id: OverlapProblemTPAnalyzer.cc,v 1.1 2011/01/22 17:59:36 venturia Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 #include <numeric>
24 #include <vector>
25 
26 // user include files
29 
34 
36 
38 
41 
43 
46 
50 
52 
53 #include "TH1F.h"
54 #include "TH2F.h"
55 //
56 // class decleration
57 //
58 
59 
61 public:
63  ~OverlapProblemTPAnalyzer() override;
64 
65 private:
66  void beginRun(const edm::Run&, const edm::EventSetup&) override;
67  void endRun(const edm::Run&, const edm::EventSetup&) override;
68  void analyze(const edm::Event&, const edm::EventSetup&) override;
69 
70  // ----------member data ---------------------------
71 
72  TH1F* m_ptp;
73  TH1F* m_etatp;
74  TH1F* m_nhits;
75  TH1F* m_nrechits;
77  TH1F* m_nassotk;
78  TH1F* m_pdgid;
79  TH1F* m_llbit;
80  TH1F* m_statustp;
81 
82  std::vector<TH1F*> m_simhitytecr;
83  std::vector<TH1F*> m_assosimhitytecr;
84 
85 
89 };
90 
91 //
92 // constants, enums and typedefs
93 //
94 
95 //
96 // static data member definitions
97 //
98 
99 //
100 // constructors and destructor
101 //
104  m_tpcollToken(consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("trackingParticlesCollection"))),
105  m_trkcollToken(consumes<edm::View<reco::Track> >(iConfig.getParameter<edm::InputTag>("trackCollection"))),
106  m_associatorToken(consumes<reco::TrackToTrackingParticleAssociator>(edm::InputTag("trackAssociatorByHits")))
107 {
108  //now do what ever initialization is needed
109 
110 
111 
113 
114  m_ptp = tfserv->make<TH1F>("tpmomentum","Tracking Particle momentum",100,0.,200.);
115  m_etatp = tfserv->make<TH1F>("tpeta","Tracking Particle pseudorapidity",100,-4.,4.);
116  m_nhits = tfserv->make<TH1F>("nhits","Tracking Particle associated hits",100,-0.5,99.5);
117  m_nrechits = tfserv->make<TH1F>("nrechits","Tracking Particle associated rec hits",100,-0.5,99.5);
118  m_nrecvssimhits = tfserv->make<TH2F>("nrecvssimhits","Tracking Particle associated rec hits vs sim hits",
119  100,-0.5,99.5,100,-0.5,99.5);
120  m_nassotk = tfserv->make<TH1F>("nassotk","Number of assocated reco tracks",10,-0.5,9.5);
121 
122  m_pdgid = tfserv->make<TH1F>("pdgid","Tracking Particle PDG id",1000,-500.5,499.5);
123  m_llbit = tfserv->make<TH1F>("llbit","Tracking Particle LongLived bit",2,-0.5,1.5);
124  m_statustp = tfserv->make<TH1F>("statustp","Tracking Particle status",2000,-1000.5,999.5);
125 
126  for(unsigned int ring=0;ring<7;++ring) {
127 
128  char name[100];
129  char title[100];
130 
131  sprintf(name,"simytecr_%d",ring+1);
132  sprintf(title,"SimHit local Y TEC ring %d",ring+1);
133 
134  m_simhitytecr.push_back(tfserv->make<TH1F>(name,title,200,-20.,20.));
135 
136  sprintf(name,"assosimytecr_%d",ring+1);
137  sprintf(title,"SimHit local Y TEC ring %d with associated RecHit",ring+1);
138 
139  m_assosimhitytecr.push_back(tfserv->make<TH1F>(name,title,200,-20.,20.));
140 
141  }
142 
143 }
144 
145 
147 {
148 
149  // do anything here that needs to be done at desctruction time
150  // (e.g. close files, deallocate resources etc.)
151 
152 }
153 
154 
155 //
156 // member functions
157 //
158 
159 // ------------ method called to for each event ------------
160 void
162 {
163  using namespace edm;
164 
165  // reco track Handle
166 
167  // Handle<reco::TrackCollection> trkcoll;
169  iEvent.getByToken(m_trkcollToken,trkcoll);
170 
172  iEvent.getByToken(m_tpcollToken,tpcoll);
173 
175  iEvent.getByToken(m_associatorToken,tahandle);
176 
177  // associate reco to sim tracks
178 
179  reco::SimToRecoCollection srcoll = tahandle->associateSimToReco(trkcoll,tpcoll);
180 
181  // loop on Handle with index and use find
182 
183  for(unsigned int index=0 ; index != tpcoll->size() ; ++index) {
184 
185  // get TrackingParticleRef
186 
187  const TrackingParticleRef tp(tpcoll,index);
188 
189  if(std::abs(tp->pdgId())!=13) continue;
190 
191  // get the SimHIt from tracker only
192 
193  // Commented since the new TP's do not have this method
194  // std::vector<PSimHit> tksimhits = tp->trackPSimHit(DetId::Tracker);
195 
196 
197  m_ptp->Fill(tp->p());
198  m_etatp->Fill(tp->eta());
199  // m_nhits->Fill(tp->matchedHit());
200  // With the new Tracking Particles I have to use a different method
201  // m_nhits->Fill(tksimhits.size());
202  m_nhits->Fill(tp->numberOfTrackerHits());
203 
204 
205  m_pdgid->Fill(tp->pdgId());
206  m_llbit->Fill(tp->longLived());
207  m_statustp->Fill(tp->status());
208 
209  // prepare a vector of TrackingRecHit from the associated reco tracks
210 
211  std::vector<DetId> rechits;
212 
213  // look at associated tracks
214 
215  if(srcoll.find(tp) != srcoll.end()) {
216  reco::SimToRecoCollection::result_type trks = srcoll[tp];
217  m_nassotk->Fill(trks.size());
218 
219  // loop on associated tracks and fill TrackingRecHit vector
220  for(reco::SimToRecoCollection::result_type::const_iterator trk = trks.begin();trk!=trks.end();++trk) {
221  for(trackingRecHit_iterator rh = trk->first->recHitsBegin() ; rh!=trk->first->recHitsEnd() ; ++rh) {
222  rechits.push_back((*rh)->geographicalId());
223  }
224 
225  }
226 
227  }
228  else {
229  m_nassotk->Fill(0.);
230  edm::LogInfo("NoAssociatedRecoTrack") << "No associated reco track for TP with p = " << tp->p()
231  << " and eta = " << tp->eta() ;
232  }
233 
234  m_nrechits->Fill(rechits.size());
235  // new method used to be compliant with the new TP's
236  m_nrecvssimhits->Fill(tp->numberOfTrackerHits(),rechits.size());
237 
238  LogDebug("RecHitDetId") << "List of " << rechits.size() << " rechits detid from muon with p = " << tp->p()
239  << "and eta = " << tp->eta();
240  for(unsigned int i=0;i<rechits.size();++i) {
241  LogTrace("RecHitDetId") << rechits[i].rawId();
242  }
243 
244 
245  // loop on sim hits
246 
247 
248  LogDebug("SimHitDetId") << "List of " << tp->numberOfTrackerHits() << " simhits detid from muon with p = " << tp->p()
249  << "and eta = " << tp->eta();
250 
251  }
252 
253 }
254 
255 void
257 {
258 }
259 
260 void
262 {
263 }
264 
265 
266 
267 //define this as a plug-in
#define LogDebug(id)
void endRun(const edm::Run &, const edm::EventSetup &) override
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< TrackingParticle > TrackingParticleCollection
const_iterator end() const
last iterator over the map (read only)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const_iterator find(const key_type &k) const
find element with specified reference key
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void beginRun(const edm::Run &, const edm::EventSetup &) override
edm::EDGetTokenT< TrackingParticleCollection > m_tpcollToken
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< edm::View< reco::Track > > m_trkcollToken
int iEvent
Definition: GenABIO.cc:230
OverlapProblemTPAnalyzer(const edm::ParameterSet &)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define LogTrace(id)
std::vector< TH1F * > m_simhitytecr
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > m_associatorToken
fixed size matrix
HLT enums.
std::vector< TH1F * > m_assosimhitytecr
reco::SimToRecoCollection associateSimToReco(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
Definition: Run.h:43