CMS 3D CMS Logo

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