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