CMS 3D CMS Logo

DTHitAssociator.cc
Go to the documentation of this file.
11 #include <sstream>
12 #include <string>
13 
14 using namespace std;
15 
16 // Constructor
18  : DTsimhitsTag(conf.getParameter<edm::InputTag>("DTsimhitsTag")),
19  DTsimhitsXFTag(conf.getParameter<edm::InputTag>("DTsimhitsXFTag")),
20  DTdigiTag(conf.getParameter<edm::InputTag>("DTdigiTag")),
21  DTdigisimlinkTag(conf.getParameter<edm::InputTag>("DTdigisimlinkTag")),
22  DTrechitTag(conf.getParameter<edm::InputTag>("DTrechitTag")),
23  geomToken(iC.esConsumes()),
24  // nice printout of DT hits
25  dumpDT(conf.getParameter<bool>("dumpDT")),
26  // CrossingFrame used or not ?
27  crossingframe(conf.getParameter<bool>("crossingframe")),
28  // Event contain the DTDigiSimLink collection ?
29  links_exist(conf.getParameter<bool>("links_exist")),
30  // associatorByWire links to a RecHit all the "valid" SimHits on the same
31  // DT wire
32  associatorByWire(conf.getParameter<bool>("associatorByWire")) {
33  if (crossingframe) {
35  } else if (!DTsimhitsTag.label().empty()) {
37  }
40 
41  if (dumpDT) {
43  }
44 
45  if (!links_exist && !associatorByWire) {
46  LogTrace("DTHitAssociator") << "*** WARNING in DTHitAssociator::DTHitAssociator: associatorByWire "
47  "reset to TRUE !"
48  << " \t (missing DTDigiSimLinkCollection ?)";
49  associatorByWire = true;
50  }
51 }
52 
54  const edm::EventSetup &iSetup,
55  const Config &conf,
56  bool printRtS)
57  : // input collection labels
58  config_(conf),
59  printRtS(true)
60 
61 {
62  initEvent(iEvent, iSetup);
63 }
64 
66  LogTrace("DTHitAssociator") << "DTHitAssociator constructor: dumpDT = " << config_.dumpDT
67  << ", crossingframe = " << config_.crossingframe
68  << ", links_exist = " << config_.links_exist
69  << ", associatorByWire = " << config_.associatorByWire;
70 
71  // need DT Geometry to discard hits for which the drift time parametrisation
72  // is not applicable
74 
75  // Get the DT SimHits from the event and map PSimHit by DTWireId
76  mapOfSimHit.clear();
77  bool takeHit(true);
78 
79  if (config_.crossingframe) {
80  LogTrace("DTHitAssociator") << "getting CrossingFrame<PSimHit> collection - " << config_.DTsimhitsXFTag;
82  unique_ptr<MixCollection<PSimHit>> DTsimhits(new MixCollection<PSimHit>(&xFrame));
83  LogTrace("DTHitAssociator") << "... size = " << DTsimhits->size();
85  for (isimhit = DTsimhits->begin(); isimhit != DTsimhits->end(); isimhit++) {
86  DTWireId wireid((*isimhit).detUnitId());
87  takeHit = SimHitOK(muonGeom, *isimhit);
88  mapOfSimHit[wireid].push_back(make_pair(*isimhit, takeHit));
89  }
90  } else if (!config_.DTsimhitsTag.label().empty()) {
91  LogTrace("DTHitAssociator") << "getting PSimHit collection - " << config_.DTsimhitsTag;
92  edm::PSimHitContainer const &DTsimhits = iEvent.get(config_.DTsimhitsToken);
93  LogTrace("DTHitAssociator") << "... size = " << DTsimhits.size();
94  edm::PSimHitContainer::const_iterator isimhit;
95  for (isimhit = DTsimhits.begin(); isimhit != DTsimhits.end(); isimhit++) {
96  DTWireId wireid((*isimhit).detUnitId());
97  takeHit = SimHitOK(muonGeom, *isimhit);
98  mapOfSimHit[wireid].push_back(make_pair(*isimhit, takeHit));
99  }
100  }
101 
102  // Get the DT Digi collection from the event
103  mapOfDigi.clear();
104  LogTrace("DTHitAssociator") << "getting DTDigi collection - " << config_.DTdigiTag;
106 
107  if (digis.isValid()) {
108  // Map DTDigi by DTWireId
109  for (DTDigiCollection::DigiRangeIterator detUnit = digis->begin(); detUnit != digis->end(); ++detUnit) {
110  const DTLayerId &layerid = (*detUnit).first;
111  const DTDigiCollection::Range &range = (*detUnit).second;
112 
114  for (digi = range.first; digi != range.second; ++digi) {
115  DTWireId wireid(layerid, (*digi).wire());
116  mapOfDigi[wireid].push_back(*digi);
117  }
118  }
119  } else {
120  LogTrace("DTHitAssociator") << "... NO DTDigi collection found !";
121  }
122 
123  mapOfLinks.clear();
124  if (config_.links_exist) {
125  // Get the DT DigiSimLink collection from the event and map DTDigiSimLink by
126  // DTWireId
127  LogTrace("DTHitAssociator") << "getting DTDigiSimLink collection - " << config_.DTdigisimlinkTag;
128  DTDigiSimLinkCollection const &digisimlinks = iEvent.get(config_.DTdigisimlinkToken);
129 
130  for (DTDigiSimLinkCollection::DigiRangeIterator detUnit = digisimlinks.begin(); detUnit != digisimlinks.end();
131  ++detUnit) {
132  const DTLayerId &layerid = (*detUnit).first;
133  const DTDigiSimLinkCollection::Range &range = (*detUnit).second;
134 
136  for (link = range.first; link != range.second; ++link) {
137  DTWireId wireid(layerid, (*link).wire());
138  mapOfLinks[wireid].push_back(*link);
139  }
140  }
141  }
142 
143  if (config_.dumpDT && printRtS) {
144  // Get the DT rechits from the event
145  LogTrace("DTHitAssociator") << "getting DTRecHit1DPair collection - " << config_.DTrechitTag;
146  DTRecHitCollection const &DTrechits = iEvent.get(config_.DTrechitToken);
147  LogTrace("DTHitAssociator") << "... size = " << DTrechits.size();
148 
149  // map DTRecHit1DPair by DTWireId
150  mapOfRecHit.clear();
152  for (rechit = DTrechits.begin(); rechit != DTrechits.end(); ++rechit) {
153  DTWireId wireid = (*rechit).wireId();
154  mapOfRecHit[wireid].push_back(*rechit);
155  }
156 
157  if (mapOfSimHit.end() != mapOfSimHit.begin()) {
158  LogTrace("DTHitAssociator") << "\n *** Dump DT PSimHit's ***";
159 
160  int jwire = 0;
161  int ihit = 0;
162 
163  for (SimHitMap::const_iterator mapIT = mapOfSimHit.begin(); mapIT != mapOfSimHit.end(); ++mapIT, jwire++) {
164  DTWireId wireid = (*mapIT).first;
165  for (vector<PSimHit_withFlag>::const_iterator hitIT = mapOfSimHit[wireid].begin();
166  hitIT != mapOfSimHit[wireid].end();
167  hitIT++, ihit++) {
168  PSimHit hit = hitIT->first;
169  LogTrace("DTHitAssociator") << "PSimHit " << ihit << ", wire " << wireid //<<", detID = "<<hit.detUnitId()
170  << ", SimTrack Id:" << hit.trackId() << "/Evt:(" << hit.eventId().event() << ","
171  << hit.eventId().bunchCrossing() << ") "
172  << ", pdg = " << hit.particleType() << ", procs = " << hit.processType();
173  }
174  }
175  } else {
176  LogTrace("DTHitAssociator") << "\n *** There are NO DT PSimHit's ***";
177  }
178 
179  if (mapOfDigi.end() != mapOfDigi.begin()) {
180  int jwire = 0;
181  int ihit = 0;
182 
183  for (DigiMap::const_iterator mapIT = mapOfDigi.begin(); mapIT != mapOfDigi.end(); ++mapIT, jwire++) {
184  LogTrace("DTHitAssociator") << "\n *** Dump DT digis ***";
185 
186  DTWireId wireid = (*mapIT).first;
187  for (vector<DTDigi>::const_iterator hitIT = mapOfDigi[wireid].begin(); hitIT != mapOfDigi[wireid].end();
188  hitIT++, ihit++) {
189  LogTrace("DTHitAssociator") << "DTDigi " << ihit << ", wire " << wireid << ", number = " << hitIT->number()
190  << ", TDC counts = " << hitIT->countsTDC();
191  }
192  }
193  } else {
194  LogTrace("DTHitAssociator") << "\n *** There are NO DTDigi's ***";
195  }
196 
197  if (mapOfRecHit.end() != mapOfRecHit.begin())
198  LogTrace("DTHitAssociator") << "\n *** Analyze DTRecHitCollection by DTWireId ***";
199 
200  int iwire = 0;
201  for (RecHitMap::const_iterator mapIT = mapOfRecHit.begin(); mapIT != mapOfRecHit.end(); ++mapIT, iwire++) {
202  DTWireId wireid = (*mapIT).first;
203  LogTrace("DTHitAssociator") << "\n======================================="
204  "============================";
205  LogTrace("DTHitAssociator") << "wire index = " << iwire << " *** DTWireId = "
206  << " (" << wireid << ")";
207 
208  if (mapOfSimHit.find(wireid) != mapOfSimHit.end()) {
209  LogTrace("DTHitAssociator") << "\n" << mapOfSimHit[wireid].size() << " SimHits (PSimHit):";
210 
211  for (vector<PSimHit_withFlag>::const_iterator hitIT = mapOfSimHit[wireid].begin();
212  hitIT != mapOfSimHit[wireid].end();
213  ++hitIT) {
214  stringstream tId;
215  tId << (hitIT->first).trackId();
216  string simhitlog = "\t SimTrack Id = " + tId.str();
217  if (hitIT->second)
218  simhitlog = simhitlog + "\t -VALID HIT-";
219  else
220  simhitlog = simhitlog + "\t -not valid hit-";
221  LogTrace("DTHitAssociator") << simhitlog;
222  }
223  }
224 
225  if (mapOfLinks.find(wireid) != mapOfLinks.end()) {
226  LogTrace("DTHitAssociator") << "\n" << mapOfLinks[wireid].size() << " Links (DTDigiSimLink):";
227 
228  for (vector<DTDigiSimLink>::const_iterator hitIT = mapOfLinks[wireid].begin();
229  hitIT != mapOfLinks[wireid].end();
230  ++hitIT) {
231  LogTrace("DTHitAssociator") << "\t digi number = " << hitIT->number() << ", time = " << hitIT->time()
232  << ", SimTrackId = " << hitIT->SimTrackId();
233  }
234  }
235 
236  if (mapOfDigi.find(wireid) != mapOfDigi.end()) {
237  LogTrace("DTHitAssociator") << "\n" << mapOfDigi[wireid].size() << " Digis (DTDigi):";
238  for (vector<DTDigi>::const_iterator hitIT = mapOfDigi[wireid].begin(); hitIT != mapOfDigi[wireid].end();
239  ++hitIT) {
240  LogTrace("DTHitAssociator") << "\t digi number = " << hitIT->number() << ", time = " << hitIT->time();
241  }
242  }
243 
244  LogTrace("DTHitAssociator") << "\n" << (*mapIT).second.size() << " RecHits (DTRecHit1DPair):";
245  for (vector<DTRecHit1DPair>::const_iterator vIT = (*mapIT).second.begin(); vIT != (*mapIT).second.end(); ++vIT) {
246  LogTrace("DTHitAssociator") << "\t digi time = " << vIT->digiTime();
247  }
248  }
249  }
250 }
251 // end of constructor
252 
253 std::vector<DTHitAssociator::SimHitIdpr> DTHitAssociator::associateHitId(const TrackingRecHit &hit) const {
254  std::vector<SimHitIdpr> simtrackids;
255  const TrackingRecHit *hitp = &hit;
256  const DTRecHit1D *dtrechit = dynamic_cast<const DTRecHit1D *>(hitp);
257 
258  if (dtrechit) {
259  simtrackids = associateDTHitId(dtrechit);
260  } else {
261  LogTrace("DTHitAssociator") << "*** WARNING in DTHitAssociator::associateHitId, null dynamic_cast "
262  "!";
263  }
264  return simtrackids;
265 }
266 
267 std::vector<DTHitAssociator::SimHitIdpr> DTHitAssociator::associateDTHitId(const DTRecHit1D *dtrechit) const {
268  std::vector<SimHitIdpr> matched;
269 
270  DTWireId wireid = dtrechit->wireId();
271 
273  // matching based on DTWireId : take only "valid" SimHits on that wire
274 
275  auto found = mapOfSimHit.find(wireid);
276  if (found != mapOfSimHit.end()) {
277  for (vector<PSimHit_withFlag>::const_iterator hitIT = found->second.begin(); hitIT != found->second.end();
278  ++hitIT) {
279  bool valid_hit = hitIT->second;
280  PSimHit this_hit = hitIT->first;
281 
282  if (valid_hit) {
283  SimHitIdpr currentId(this_hit.trackId(), this_hit.eventId());
284  matched.push_back(currentId);
285  }
286  }
287  }
288  }
289 
290  else {
291  // matching based on DTDigiSimLink
292 
293  float theTime = dtrechit->digiTime();
294  int theNumber(-1);
295 
296  auto found = mapOfLinks.find(wireid);
297 
298  if (found != mapOfLinks.end()) {
299  // DTDigiSimLink::time() is set equal to DTDigi::time() only for the
300  // DTDigiSimLink corresponding to the digitizied PSimHit other
301  // DTDigiSimLinks associated to the same DTDigi are identified by having
302  // the same DTDigiSimLink::number()
303 
304  // first find the associated digi Number
305  for (vector<DTDigiSimLink>::const_iterator linkIT = found->second.begin(); linkIT != found->second.end();
306  ++linkIT) {
307  float digitime = linkIT->time();
308  if (fabs(digitime - theTime) < 0.1) {
309  theNumber = linkIT->number();
310  }
311  }
312 
313  // then get all the DTDigiSimLinks with that digi Number (corresponding to
314  // valid SimHits
315  // within a time window of the order of the time resolution, specified in
316  // the DTDigitizer)
317  for (vector<DTDigiSimLink>::const_iterator linkIT = found->second.begin(); linkIT != found->second.end();
318  ++linkIT) {
319  int digiNr = linkIT->number();
320  if (digiNr == theNumber) {
321  SimHitIdpr currentId(linkIT->SimTrackId(), linkIT->eventId());
322  matched.push_back(currentId);
323  }
324  }
325 
326  } else {
327  LogTrace("DTHitAssociator") << "*** ERROR in DTHitAssociator::associateDTHitId - DTRecHit1D: " << *dtrechit
328  << " has no associated DTDigiSimLink !" << endl;
329  return matched;
330  }
331  }
332 
333  return matched;
334 }
335 
336 std::vector<PSimHit> DTHitAssociator::associateHit(const TrackingRecHit &hit) const {
337  std::vector<PSimHit> simhits;
338  std::vector<SimHitIdpr> simtrackids;
339 
340  const TrackingRecHit *hitp = &hit;
341  const DTRecHit1D *dtrechit = dynamic_cast<const DTRecHit1D *>(hitp);
342 
343  if (dtrechit) {
344  DTWireId wireid = dtrechit->wireId();
345 
347  // matching based on DTWireId : take only "valid" SimHits on that wire
348 
349  auto found = mapOfSimHit.find(wireid);
350  if (found != mapOfSimHit.end()) {
351  for (vector<PSimHit_withFlag>::const_iterator hitIT = found->second.begin(); hitIT != found->second.end();
352  ++hitIT) {
353  bool valid_hit = hitIT->second;
354  if (valid_hit)
355  simhits.push_back(hitIT->first);
356  }
357  }
358  } else {
359  // matching based on DTDigiSimLink
360 
361  simtrackids = associateDTHitId(dtrechit);
362 
363  for (vector<SimHitIdpr>::const_iterator idIT = simtrackids.begin(); idIT != simtrackids.end(); ++idIT) {
364  uint32_t trId = idIT->first;
365  EncodedEventId evId = idIT->second;
366 
367  auto found = mapOfSimHit.find(wireid);
368  if (found != mapOfSimHit.end()) {
369  for (vector<PSimHit_withFlag>::const_iterator hitIT = found->second.begin(); hitIT != found->second.end();
370  ++hitIT) {
371  if (hitIT->first.trackId() == trId && hitIT->first.eventId() == evId)
372  simhits.push_back(hitIT->first);
373  }
374  }
375  }
376  }
377 
378  } else {
379  LogTrace("DTHitAssociator") << "*** WARNING in DTHitAssociator::associateHit, null dynamic_cast !";
380  }
381  return simhits;
382 }
383 
384 bool DTHitAssociator::SimHitOK(const edm::ESHandle<DTGeometry> &muongeom, const PSimHit &simhit) {
385  bool result = true;
386 
387  DTWireId wireid(simhit.detUnitId());
388  const DTLayer *dtlayer = muongeom->layer(wireid.layerId());
389  LocalPoint entryP = simhit.entryPoint();
390  LocalPoint exitP = simhit.exitPoint();
391  const DTTopology &topo = dtlayer->specificTopology();
392  float xwire = topo.wirePosition(wireid.wire());
393  float xEntry = entryP.x() - xwire;
394  float xExit = exitP.x() - xwire;
395  DTTopology::Side entrySide = topo.onWhichBorder(xEntry, entryP.y(), entryP.z());
396  DTTopology::Side exitSide = topo.onWhichBorder(xExit, exitP.y(), exitP.z());
397 
398  bool noParametrisation =
399  ((entrySide == DTTopology::none || exitSide == DTTopology::none) || (entrySide == exitSide) ||
400  ((entrySide == DTTopology::xMin && exitSide == DTTopology::xMax) ||
401  (entrySide == DTTopology::xMax && exitSide == DTTopology::xMin)));
402 
403  // discard hits where parametrization can not be used
404  if (noParametrisation) {
405  result = false;
406  return result;
407  }
408 
409  float x;
410  LocalPoint hitPos = simhit.localPosition();
411 
412  if (fabs(hitPos.z()) < 0.002) {
413  // hit center within 20 um from z=0, no need to extrapolate.
414  x = hitPos.x() - xwire;
415  } else {
416  x = xEntry - (entryP.z() * (xExit - xEntry)) / (exitP.z() - entryP.z());
417  }
418 
419  // discard hits where x is out of range of the parametrization (|x|>2.1 cm)
420  x *= 10.; // cm -> mm
421  if (fabs(x) > 21.)
422  result = false;
423 
424  return result;
425 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::InputTag DTdigisimlinkTag
std::vector< PSimHit > associateHit(const TrackingRecHit &hit) const
std::pair< uint32_t, EncodedEventId > SimHitIdpr
Config(const edm::ParameterSet &, edm::ConsumesCollector iC)
unsigned int detUnitId() const
Definition: PSimHit.h:99
T z() const
Definition: PV3DBase.h:61
edm::EDGetTokenT< DTDigiCollection > DTdigiToken
Side onWhichBorder(float x, float y, float z) const
Definition: DTTopology.cc:82
RecHitMap mapOfRecHit
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:48
std::vector< SimHitIdpr > associateDTHitId(const DTRecHit1D *dtrechit) const
std::string const & label() const
Definition: InputTag.h:36
#define LogTrace(id)
DTHitAssociator(const edm::Event &, const edm::EventSetup &, const Config &, bool printRtS)
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
SimHitMap mapOfSimHit
edm::InputTag DTsimhitsXFTag
edm::EDGetTokenT< DTRecHitCollection > DTrechitToken
edm::EDGetTokenT< DTDigiSimLinkCollection > DTdigisimlinkToken
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:45
edm::ESGetToken< DTGeometry, MuonGeometryRecord > geomToken
unsigned int trackId() const
Definition: PSimHit.h:108
EncodedEventId eventId() const
Definition: PSimHit.h:119
Side
Sides of the cell.
Definition: DTTopology.h:89
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
Config const & config_
DigiRangeIterator begin() const
edm::EDGetTokenT< edm::PSimHitContainer > DTsimhitsToken
Local3DPoint localPosition() const
Definition: PSimHit.h:54
std::pair< const_iterator, const_iterator > Range
std::vector< DigiType >::const_iterator const_iterator
bool isValid() const
Definition: HandleBase.h:70
edm::InputTag DTsimhitsTag
HLT enums.
void initEvent(const edm::Event &, const edm::EventSetup &)
Definition: Config.py:1
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &hit) const
std::vector< PSimHit > PSimHitContainer
edm::EDGetTokenT< CrossingFrame< PSimHit > > DTsimhitsXFToken
bool SimHitOK(const edm::ESHandle< DTGeometry > &, const PSimHit &)
DigiRangeIterator end() const
float wirePosition(int wireNumber) const
Returns the x position in the layer of a given wire number.
Definition: DTTopology.cc:63
const DTLayer * layer(const DTLayerId &id) const
Return a layer given its id.
Definition: DTGeometry.cc:96
A container for a generic type of digis indexed by some index, implemented with a map<IndexType...