CMS 3D CMS Logo

SiPixelPhase1HitsV.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1HitsV
4 // Class: SiPixelPhase1HitsV
5 //
6 
7 // Original Author: Marcel Schneider
8 // Additional Authors: Alexander Morton - modifying code for validation use
9 
13 
16 
20 
23 
24 //class TrackAssociatorByHits;
25 
27  SiPixelPhase1Base(iConfig),
28  pixelBarrelLowToken_ ( consumes<edm::PSimHitContainer>(iConfig.getParameter<edm::InputTag>("pixBarrelLowSrc")) ),
29  pixelBarrelHighToken_ ( consumes<edm::PSimHitContainer>(iConfig.getParameter<edm::InputTag>("pixBarrelHighSrc")) ),
30  pixelForwardLowToken_ ( consumes<edm::PSimHitContainer>(iConfig.getParameter<edm::InputTag>("pixForwardLowSrc")) ),
31  pixelForwardHighToken_ ( consumes<edm::PSimHitContainer>(iConfig.getParameter<edm::InputTag>("pixForwardHighSrc")) ),
32 
33  tracksToken_ ( consumes< edm::View<reco::Track> >(iConfig.getParameter<edm::InputTag>("tracksTag")) ),
34  tpToken_ ( consumes< TrackingParticleCollection >(iConfig.getParameter<edm::InputTag>("tpTag")) ),
35  trackAssociatorByHitsToken_ ( consumes< reco::TrackToTrackingParticleAssociator >(iConfig.getParameter<edm::InputTag>("trackAssociatorByHitsTag")) )
36 {}
37 
39 
41  iEvent.getByToken(pixelBarrelLowToken_, barrelLowInput);
42  if (!barrelLowInput.isValid()) return;
43 
44  edm::Handle<edm::PSimHitContainer> barrelHighInput;
45  iEvent.getByToken(pixelBarrelHighToken_, barrelHighInput);
46  if (!barrelHighInput.isValid()) return;
47 
48  edm::Handle<edm::PSimHitContainer> forwardLowInput;
49  iEvent.getByToken(pixelForwardLowToken_, forwardLowInput);
50  if (!forwardLowInput.isValid()) return;
51 
52  edm::Handle<edm::PSimHitContainer> forwardHighInput;
53  iEvent.getByToken(pixelForwardHighToken_, forwardHighInput);
54  if (!forwardHighInput.isValid()) return;
55 
56  edm::PSimHitContainer::const_iterator it;
57 
58  // Get geometry information
59 
61  iSetup.get<TrackerDigiGeometryRecord>().get( tracker );
62 
63  // get low barrel info
64  for (it = barrelLowInput->begin(); it != barrelLowInput->end(); ++it) {
65  auto id = DetId(it->detUnitId());
66  const GeomDetUnit * det=(const GeomDetUnit*)tracker->idToDetUnit( id );
67  GlobalPoint gpos=det->toGlobal(it->localPosition());
68 
69  float tof = it->timeOfFlight();
70  float globalR = gpos.mag();
71 
72  float energyLoss = it->energyLoss();
73 
74  float entryExitX = ( it->entryPoint().x() - it->exitPoint().x() );
75  float entryExitY = ( it->entryPoint().y() - it->exitPoint().y() );
76  float entryExitZ = std::abs( it->entryPoint().z() - it->exitPoint().z() );
77 
78  float localX = it->localPosition().x();
79  float localY = it->localPosition().y();
80  float localZ = it->localPosition().z();
81  float localPhi = it->localPosition().phi();
82  float localEta = it->localPosition().eta();
83 
84  histo[TOF_R].fill(globalR, tof, id, &iEvent);
85  histo[ELOSS].fill(energyLoss, id, &iEvent);
86  histo[ENTRY_EXIT_X].fill(entryExitX, id, &iEvent);
87  histo[ENTRY_EXIT_Y].fill(entryExitY, id, &iEvent);
88  histo[ENTRY_EXIT_Z].fill(entryExitZ, id, &iEvent);
89  histo[LOCAL_X].fill(localX, id, &iEvent);
90  histo[LOCAL_Y].fill(localY, id, &iEvent);
91  histo[LOCAL_Z].fill(localZ, id, &iEvent);
92  histo[LOCAL_PHI].fill(localPhi, id, &iEvent);
93  histo[LOCAL_ETA].fill(localEta, id, &iEvent);
94  }
95  // get high barrel info
96  for (it = barrelHighInput->begin(); it != barrelHighInput->end(); ++it) {
97  auto id = DetId(it->detUnitId());
98  const GeomDetUnit * det=(const GeomDetUnit*)tracker->idToDetUnit( id );
99  GlobalPoint gpos=det->toGlobal(it->localPosition());
100 
101  float tof = it->timeOfFlight();
102  float globalR = gpos.mag();
103 
104  float energyLoss = it->energyLoss();
105 
106  float entryExitX = ( it->entryPoint().x() - it->exitPoint().x() );
107  float entryExitY = ( it->entryPoint().y() - it->exitPoint().y() );
108  float entryExitZ = std::abs( it->entryPoint().z() - it->exitPoint().z() );
109 
110  float localX = it->localPosition().x();
111  float localY = it->localPosition().y();
112  float localZ = it->localPosition().z();
113  float localPhi = it->localPosition().phi();
114  float localEta = it->localPosition().eta();
115 
116  histo[TOF_R].fill(globalR, tof, id, &iEvent);
117  histo[ELOSS].fill(energyLoss, id, &iEvent);
118  histo[ENTRY_EXIT_X].fill(entryExitX, id, &iEvent);
119  histo[ENTRY_EXIT_Y].fill(entryExitY, id, &iEvent);
120  histo[ENTRY_EXIT_Z].fill(entryExitZ, id, &iEvent);
121  histo[LOCAL_X].fill(localX, id, &iEvent);
122  histo[LOCAL_Y].fill(localY, id, &iEvent);
123  histo[LOCAL_Z].fill(localZ, id, &iEvent);
124  histo[LOCAL_PHI].fill(localPhi, id, &iEvent);
125  histo[LOCAL_ETA].fill(localEta, id, &iEvent);
126  }
127 
128  // get low forward info
129  for (it = forwardLowInput->begin(); it != forwardLowInput->end(); ++it) {
130  auto id = DetId(it->detUnitId());
131  const GeomDetUnit * det=(const GeomDetUnit*)tracker->idToDetUnit( id );
132  GlobalPoint gpos=det->toGlobal(it->localPosition());
133 
134  float tof = it->timeOfFlight();
135  float globalR = gpos.mag();
136 
137  float energyLoss = it->energyLoss();
138 
139  float entryExitX = ( it->entryPoint().x() - it->exitPoint().x() );
140  float entryExitY = ( it->entryPoint().y() - it->exitPoint().y() );
141  float entryExitZ = std::abs( it->entryPoint().z() - it->exitPoint().z() );
142 
143  float localX = it->localPosition().x();
144  float localY = it->localPosition().y();
145  float localZ = it->localPosition().z();
146  float localPhi = it->localPosition().phi();
147  float localEta = it->localPosition().eta();
148 
149  histo[TOF_R].fill(globalR, tof, id, &iEvent);
150  histo[ELOSS].fill(energyLoss, id, &iEvent);
151  histo[ENTRY_EXIT_X].fill(entryExitX, id, &iEvent);
152  histo[ENTRY_EXIT_Y].fill(entryExitY, id, &iEvent);
153  histo[ENTRY_EXIT_Z].fill(entryExitZ, id, &iEvent);
154  histo[LOCAL_X].fill(localX, id, &iEvent);
155  histo[LOCAL_Y].fill(localY, id, &iEvent);
156  histo[LOCAL_Z].fill(localZ, id, &iEvent);
157  histo[LOCAL_PHI].fill(localPhi, id, &iEvent);
158  histo[LOCAL_ETA].fill(localEta, id, &iEvent);
159  }
160 
161  // get high forward info
162  for (it = forwardHighInput->begin(); it != forwardHighInput->end(); ++it) {
163  auto id = DetId(it->detUnitId());
164  const GeomDetUnit * det=(const GeomDetUnit*)tracker->idToDetUnit( id );
165  GlobalPoint gpos=det->toGlobal(it->localPosition());
166 
167  float tof = it->timeOfFlight();
168  float globalR = gpos.mag();
169 
170  float energyLoss = it->energyLoss();
171 
172  float entryExitX = ( it->entryPoint().x() - it->exitPoint().x() );
173  float entryExitY = ( it->entryPoint().y() - it->exitPoint().y() );
174  float entryExitZ = std::abs( it->entryPoint().z() - it->exitPoint().z() );
175 
176  float localX = it->localPosition().x();
177  float localY = it->localPosition().y();
178  float localZ = it->localPosition().z();
179  float localPhi = it->localPosition().phi();
180  float localEta = it->localPosition().eta();
181 
182  histo[TOF_R].fill(globalR, tof, id, &iEvent);
183  histo[ELOSS].fill(energyLoss, id, &iEvent);
184  histo[ENTRY_EXIT_X].fill(entryExitX, id, &iEvent);
185  histo[ENTRY_EXIT_Y].fill(entryExitY, id, &iEvent);
186  histo[ENTRY_EXIT_Z].fill(entryExitZ, id, &iEvent);
187  histo[LOCAL_X].fill(localX, id, &iEvent);
188  histo[LOCAL_Y].fill(localY, id, &iEvent);
189  histo[LOCAL_Z].fill(localZ, id, &iEvent);
190  histo[LOCAL_PHI].fill(localPhi, id, &iEvent);
191  histo[LOCAL_ETA].fill(localEta, id, &iEvent);
192  }
193 
194 
195  // Sim Hit efficiency info
196  edm::Handle< edm::View<reco::Track> > trackCollectionH;
197  iEvent.getByToken(tracksToken_, trackCollectionH);
198  const edm::View<reco::Track>& tC = *(trackCollectionH.product());
199 
201  iEvent.getByToken(tpToken_,TPCollectionH);
202 
204  iEvent.getByToken(trackAssociatorByHitsToken_,theHitsAssociator);
205  if ( ! theHitsAssociator.isValid() ) {
206  throw cms::Exception ("NO VALID HIT ASSOCIATOR");
207  }
208  reco::TrackToTrackingParticleAssociator const *associatorByHits = theHitsAssociator.product();
209 
210  if ( TPCollectionH.isValid() && trackCollectionH.isValid() ) {
211  reco::RecoToSimCollection const& p = associatorByHits->associateRecoToSim(trackCollectionH,TPCollectionH);
212 
213  for(edm::View<reco::Track>::size_type i=0; i<tC.size(); ++i) {
214  edm::RefToBase<reco::Track> track(trackCollectionH, i);
215 // const reco::Track& t = *track;
216  auto id = DetId(track->innerDetId()); // histo manager requires a det ID, use innermost ID for ease
217 
218  auto iter = p.find(track);
219  histo[EFFICIENCY_TRACK].fill(iter != p.end()? 1: 0, id, &iEvent);
220  }
221 
222  }
223 
224 }
225 
227 
unsigned int size_type
Definition: View.h:90
std::vector< TrackingParticle > TrackingParticleCollection
edm::EDGetTokenT< TrackingParticleCollection > tpToken_
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
edm::EDGetTokenT< edm::PSimHitContainer > pixelBarrelHighToken_
size_type size() const
edm::EDGetTokenT< edm::PSimHitContainer > pixelForwardHighToken_
SiPixelPhase1HitsV(const edm::ParameterSet &conf)
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
int iEvent
Definition: GenABIO.cc:230
T mag() const
Definition: PV3DBase.h:67
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > trackAssociatorByHitsToken_
edm::EDGetTokenT< edm::PSimHitContainer > pixelBarrelLowToken_
bool isValid() const
Definition: HandleBase.h:74
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
reco::RecoToSimCollection associateRecoToSim(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
void analyze(const edm::Event &, const edm::EventSetup &) override
const T & get() const
Definition: EventSetup.h:59
std::vector< HistogramManager > histo
fixed size matrix
HLT enums.
std::vector< PSimHit > PSimHitContainer
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:99
edm::EDGetTokenT< edm::View< reco::Track > > tracksToken_
edm::EDGetTokenT< edm::PSimHitContainer > pixelForwardLowToken_