CMS 3D CMS Logo

DTSimHitMatcher.cc
Go to the documentation of this file.
2 
3 using namespace std;
4 
7  : MuonSimHitMatcher(ps, std::move(iC)) {
8  simHitPSet_ = ps.getParameterSet("dtSimHit");
9  verbose_ = simHitPSet_.getParameter<int>("verbose");
10  simMuOnly_ = simHitPSet_.getParameter<bool>("simMuOnly");
11  discardEleHits_ = simHitPSet_.getParameter<bool>("discardEleHits");
12 
13  simHitInput_ = iC.consumes<edm::PSimHitContainer>(
15 }
16 
19  const edm::EventSetup& iSetup) {
20  iSetup.get<MuonGeometryRecord>().get(dt_geom_);
21  if (dt_geom_.isValid()) {
22  geometry_ = &*dt_geom_;
23  } else {
24  hasGeometry_ = false;
25  edm::LogWarning("DTSimHitMatcher")
26  << "+++ Info: DT geometry is unavailable. +++\n";
27  }
28  MuonSimHitMatcher::init(iEvent, iSetup);
29 }
30 
32 void DTSimHitMatcher::match(const SimTrack& track, const SimVertex& vertex) {
33  // instantiates the track ids and simhits
34  MuonSimHitMatcher::match(track, vertex);
35 
36  if (hasGeometry_) {
38 
39  if (verbose_) {
40  edm::LogInfo("DTSimHitMatcher")
41  << "nTrackIds " << track_ids_.size() << " nSelectedDTSimHits "
42  << hits_.size() << endl;
43  edm::LogInfo("DTSimHitMatcher")
44  << "detids DT " << detIds(0).size() << endl;
45 
46  const auto& dt_det_ids = detIds(0);
47  for (const auto& id : dt_det_ids) {
48  const auto& dt_simhits = MuonSimHitMatcher::hitsInDetId(id);
49  const auto& dt_simhits_gp = simHitsMeanPosition(dt_simhits);
50  edm::LogInfo("DTSimHitMatcher")
51  << "DTWireId " << DTWireId(id) << ": nHits " << dt_simhits.size()
52  << " eta " << dt_simhits_gp.eta() << " phi " << dt_simhits_gp.phi()
53  << " nCh " << chamber_to_hits_[id].size() << endl;
54  }
55  }
56  }
57 }
58 
60  for (const auto& track_id : track_ids_) {
61  for (const auto& h : simHits_) {
62  if (h.trackId() != track_id) continue;
63  int pdgid = h.particleType();
64  if (simMuOnly_ && std::abs(pdgid) != 13) continue;
65  // discard electron hits in the DT chambers
66  if (discardEleHits_ && pdgid == 11) continue;
67 
68  const DTWireId layer_id(h.detUnitId());
69  detid_to_hits_[h.detUnitId()].push_back(h);
70  hits_.push_back(h);
71  layer_to_hits_[layer_id.layerId().rawId()].push_back(h);
72  superlayer_to_hits_[layer_id.superlayerId().rawId()].push_back(h);
73  chamber_to_hits_[layer_id.chamberId().rawId()].push_back(h);
74  }
75  }
76 }
77 
78 std::set<unsigned int> DTSimHitMatcher::detIds(int dt_type) const {
79  std::set<unsigned int> result;
80  for (const auto& p : detid_to_hits_) {
81  const auto& id = p.first;
82  if (dt_type > 0) {
83  DTWireId detId(id);
84  if (MuonHitHelper::toDTType(detId.wheel(), detId.station()) != dt_type)
85  continue;
86  }
87  result.insert(id);
88  }
89  return result;
90 }
91 
92 std::set<unsigned int> DTSimHitMatcher::chamberIds(int dt_type) const {
93  std::set<unsigned int> result;
94  for (const auto& p : chamber_to_hits_) {
95  const auto& id = p.first;
96  if (dt_type > 0) {
97  DTChamberId detId(id);
98  if (MuonHitHelper::toDTType(detId.wheel(), detId.station()) != dt_type)
99  continue;
100  }
101  result.insert(id);
102  }
103  return result;
104 }
105 
106 std::set<unsigned int> DTSimHitMatcher::layerIds() const {
107  std::set<unsigned int> result;
108  for (const auto& p : layer_to_hits_) result.insert(p.first);
109  return result;
110 }
111 
112 std::set<unsigned int> DTSimHitMatcher::superlayerIds() const {
113  std::set<unsigned int> result;
114  for (const auto& p : superlayer_to_hits_) result.insert(p.first);
115  return result;
116 }
117 
119  unsigned int detid) const {
120  if (!MuonHitHelper::isDT(detid)) return no_hits_;
121 
122  const DTWireId id(detid);
123  if (layer_to_hits_.find(id.layerId().rawId()) == layer_to_hits_.end())
124  return no_hits_;
125  return layer_to_hits_.at(id.layerId().rawId());
126 }
127 
129  unsigned int detid) const {
130  if (!MuonHitHelper::isDT(detid)) return no_hits_;
131 
132  const DTWireId id(detid);
133  if (superlayer_to_hits_.find(id.superlayerId().rawId()) ==
134  superlayer_to_hits_.end())
135  return no_hits_;
136  return superlayer_to_hits_.at(id.superlayerId().rawId());
137 }
138 
140  unsigned int detid) const {
141  if (!MuonHitHelper::isDT(detid)) return no_hits_;
142 
143  const DTWireId id(detid);
144  if (chamber_to_hits_.find(id.chamberId().rawId()) == chamber_to_hits_.end())
145  return no_hits_;
146  return chamber_to_hits_.at(id.chamberId().rawId());
147 }
148 
149 bool DTSimHitMatcher::hitStation(int st, int nsuperlayers, int nlayers) const {
150  int nst = 0;
151  for (const auto& ddt : chamberIds()) {
152  const DTChamberId id(ddt);
153  if (id.station() != st) continue;
154 
155  // require at least 1 superlayer
156  const int nsl(nSuperLayersWithHitsInChamber(id.rawId()));
157  if (nsl < nsuperlayers) continue;
158 
159  // require at least 3 layers hit per chamber
160  const int nl(nLayersWithHitsInChamber(id.rawId()));
161  if (nl < nlayers) continue;
162  ++nst;
163  }
164  return nst;
165 }
166 
167 int DTSimHitMatcher::nStations(int nsuperlayers, int nlayers) const {
168  return (hitStation(1, nsuperlayers, nlayers) +
169  hitStation(2, nsuperlayers, nlayers) +
170  hitStation(3, nsuperlayers, nlayers) +
171  hitStation(4, nsuperlayers, nlayers));
172 }
173 
174 int DTSimHitMatcher::nCellsWithHitsInLayer(unsigned int detid) const {
175  set<int> layers_with_hits;
176  const auto& hits = hitsInLayer(detid);
177  for (const auto& h : hits) {
178  if (MuonHitHelper::isDT(detid)) {
179  const DTWireId idd(h.detUnitId());
180  layers_with_hits.insert(idd.wire());
181  }
182  }
183  return layers_with_hits.size();
184 }
185 
186 int DTSimHitMatcher::nLayersWithHitsInSuperLayer(unsigned int detid) const {
187  set<int> layers_with_hits;
188  const auto& hits = hitsInSuperLayer(detid);
189  for (const auto& h : hits) {
190  if (MuonHitHelper::isDT(detid)) {
191  const DTLayerId idd(h.detUnitId());
192  layers_with_hits.insert(idd.layer());
193  }
194  }
195  return layers_with_hits.size();
196 }
197 
198 int DTSimHitMatcher::nSuperLayersWithHitsInChamber(unsigned int detid) const {
199  set<int> sl_with_hits;
200  const auto& hits = MuonSimHitMatcher::hitsInChamber(detid);
201  for (const auto& h : hits) {
202  if (MuonHitHelper::isDT(detid)) {
203  const DTSuperLayerId idd(h.detUnitId());
204  sl_with_hits.insert(idd.superLayer());
205  }
206  }
207  return sl_with_hits.size();
208 }
209 
210 int DTSimHitMatcher::nLayersWithHitsInChamber(unsigned int detid) const {
211  int nLayers = 0;
212  const auto& superLayers(dynamic_cast<const DTGeometry*>(geometry_)
213  ->chamber(DTChamberId(detid))
214  ->superLayers());
215  for (const auto& sl : superLayers) {
216  nLayers += nLayersWithHitsInSuperLayer(sl->id().rawId());
217  }
218  return nLayers;
219 }
221  const edm::PSimHitContainer& sim_hits) const {
222  if (sim_hits.empty()) return -1.f;
223 
224  float sums = 0.f;
225  size_t n = 0;
226  for (const auto& h : sim_hits) {
227  const LocalPoint& lp = h.entryPoint();
228  float s;
229  const auto& d = h.detUnitId();
230  if (MuonHitHelper::isDT(d)) {
231  // find nearest wire
232  s = dynamic_cast<const DTGeometry*>(geometry_)
233  ->layer(DTLayerId(d))
234  ->specificTopology()
235  .channel(lp);
236  } else
237  continue;
238  sums += s;
239  ++n;
240  }
241  if (n == 0) return -1.f;
242  return sums / n;
243 }
244 
246  unsigned int detid, int margin_n_wires) const {
247  set<unsigned int> result;
248  if (MuonHitHelper::isDT(detid)) {
249  DTLayerId id(detid);
250  int max_nwires = dynamic_cast<const DTGeometry*>(geometry_)
251  ->layer(id)
252  ->specificTopology()
253  .channels();
254  for (int wn = 0; wn <= max_nwires; ++wn) {
255  DTWireId wid(id, wn);
256  for (const auto& h : MuonSimHitMatcher::hitsInDetId(wid.rawId())) {
257  if (verbose_)
258  edm::LogInfo("DTSimHitMatcher")
259  << "central DTWireId " << wid << " simhit " << h << endl;
260  int smin = wn - margin_n_wires;
261  smin = (smin > 0) ? smin : 1;
262  int smax = wn + margin_n_wires;
263  smax = (smax <= max_nwires) ? smax : max_nwires;
264  for (int ss = smin; ss <= smax; ++ss) {
265  DTWireId widd(id, ss);
266  if (verbose_)
267  edm::LogInfo("DTSimHitMatcher")
268  << "\tadding DTWireId to collection " << widd << endl;
269  result.insert(widd.rawId());
270  }
271  }
272  }
273  }
274  return result;
275 }
276 
278  unsigned int detid, int margin_n_wires) const {
279  set<unsigned int> result;
280  const auto& layers(dynamic_cast<const DTGeometry*>(geometry_)
281  ->superLayer(DTSuperLayerId(detid))
282  ->layers());
283  for (const auto& l : layers) {
284  if (verbose_)
285  edm::LogInfo("DTSimHitMatcher")
286  << "hitWiresInDTSuperLayerId::l id " << l->id() << endl;
287  const auto& p(hitWiresInDTLayerId(l->id().rawId(), margin_n_wires));
288  result.insert(p.begin(), p.end());
289  }
290  return result;
291 }
292 
294  unsigned int detid, int margin_n_wires) const {
295  set<unsigned int> result;
296  const auto& superLayers(dynamic_cast<const DTGeometry*>(geometry_)
297  ->chamber(DTChamberId(detid))
298  ->superLayers());
299  for (const auto& sl : superLayers) {
300  if (verbose_)
301  edm::LogInfo("DTSimHitMatcher")
302  << "hitWiresInDTChamberId::sl id " << sl->id() << endl;
303  const auto& p(hitWiresInDTSuperLayerId(sl->id().rawId(), margin_n_wires));
304  result.insert(p.begin(), p.end());
305  }
306  return result;
307 }
308 
310  const std::set<unsigned int>& set) const {
311  for (const auto& p : set) {
312  DTChamberId detId(p);
313  edm::LogInfo("DTSimHitMatcher") << " " << detId << "\n";
314  }
315 }
316 
317 std::set<unsigned int> DTSimHitMatcher::chamberIdsStation(int station) const {
318  set<unsigned int> result;
319  switch (station) {
320  case 1: {
321  const auto& p1(chamberIds(MuonHitHelper::DT_MB21p));
322  const auto& p2(chamberIds(MuonHitHelper::DT_MB11p));
323  const auto& p3(chamberIds(MuonHitHelper::DT_MB01));
324  const auto& p4(chamberIds(MuonHitHelper::DT_MB11n));
325  const auto& p5(chamberIds(MuonHitHelper::DT_MB21n));
326  result.insert(p1.begin(), p1.end());
327  result.insert(p2.begin(), p2.end());
328  result.insert(p3.begin(), p3.end());
329  result.insert(p4.begin(), p4.end());
330  result.insert(p5.begin(), p5.end());
331  break;
332  }
333  case 2: {
334  const auto& p1(chamberIds(MuonHitHelper::DT_MB22p));
335  const auto& p2(chamberIds(MuonHitHelper::DT_MB12p));
336  const auto& p3(chamberIds(MuonHitHelper::DT_MB02));
337  const auto& p4(chamberIds(MuonHitHelper::DT_MB12n));
338  const auto& p5(chamberIds(MuonHitHelper::DT_MB22n));
339  result.insert(p1.begin(), p1.end());
340  result.insert(p2.begin(), p2.end());
341  result.insert(p3.begin(), p3.end());
342  result.insert(p4.begin(), p4.end());
343  result.insert(p5.begin(), p5.end());
344  break;
345  }
346  case 3: {
347  const auto& p1(chamberIds(MuonHitHelper::DT_MB23p));
348  const auto& p2(chamberIds(MuonHitHelper::DT_MB13p));
349  const auto& p3(chamberIds(MuonHitHelper::DT_MB03));
350  const auto& p4(chamberIds(MuonHitHelper::DT_MB13n));
351  const auto& p5(chamberIds(MuonHitHelper::DT_MB23n));
352  result.insert(p1.begin(), p1.end());
353  result.insert(p2.begin(), p2.end());
354  result.insert(p3.begin(), p3.end());
355  result.insert(p4.begin(), p4.end());
356  result.insert(p5.begin(), p5.end());
357  break;
358  }
359  case 4: {
360  const auto& p1(chamberIds(MuonHitHelper::DT_MB24p));
361  const auto& p2(chamberIds(MuonHitHelper::DT_MB14p));
362  const auto& p3(chamberIds(MuonHitHelper::DT_MB04));
363  const auto& p4(chamberIds(MuonHitHelper::DT_MB14n));
364  const auto& p5(chamberIds(MuonHitHelper::DT_MB24n));
365  result.insert(p1.begin(), p1.end());
366  result.insert(p2.begin(), p2.end());
367  result.insert(p3.begin(), p3.end());
368  result.insert(p4.begin(), p4.end());
369  result.insert(p5.begin(), p5.end());
370  break;
371  }
372  };
373  return result;
374 }
static bool isDT(unsigned int detId)
check detid type
Definition: MuonHitHelper.cc:3
T getParameter(std::string const &) const
const TrackingGeometry * geometry_
std::set< unsigned int > detIds(int type=MuonHitHelper::DT_ALL) const
edm::PSimHitContainer hits_
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
std::map< unsigned int, edm::PSimHitContainer > detid_to_hits_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
const edm::PSimHitContainer & hitsInChamber(unsigned int) const
void matchSimHitsToSimTrack()
int nLayersWithHitsInSuperLayer(unsigned int) const
std::set< unsigned int > hitWiresInDTSuperLayerId(unsigned int, int margin_n_wires=0) const
std::set< unsigned int > chamberIdsStation(int station) const
std::set< unsigned int > hitWiresInDTChamberId(unsigned int, int margin_n_wires=0) const
std::vector< unsigned > track_ids_
DTSimHitMatcher(const edm::ParameterSet &iPS, edm::ConsumesCollector &&iC)
std::set< unsigned int > chamberIds(int type=MuonHitHelper::DT_ALL) const
int iEvent
Definition: GenABIO.cc:224
GlobalPoint simHitsMeanPosition(const edm::PSimHitContainer &sim_hits) const
edm::PSimHitContainer simHits_
const edm::PSimHitContainer & hitsInSuperLayer(unsigned int) const
edm::ParameterSet simHitPSet_
double p4[4]
Definition: TauolaWrapper.h:92
void match(const SimTrack &t, const SimVertex &v)
do the matching
std::map< unsigned int, edm::PSimHitContainer > chamber_to_hits_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int nStations(int nsl=1, int nl=3) const
const edm::PSimHitContainer & hitsInLayer(unsigned int) const
double f[11][100]
double p2[4]
Definition: TauolaWrapper.h:90
edm::ESHandle< DTGeometry > dt_geom_
std::set< unsigned int > hitWiresInDTLayerId(unsigned int, int margin_n_wires=0) const
edm::PSimHitContainer no_hits_
int nLayersWithHitsInChamber(unsigned int) const
std::set< unsigned int > layerIds() const
bool hitStation(int, int, int) const
ParameterSet const & getParameterSet(std::string const &) const
const edm::PSimHitContainer & hitsInDetId(unsigned int) const
void match(const SimTrack &t, const SimVertex &v)
do the matching
edm::EDGetTokenT< edm::PSimHitContainer > simHitInput_
std::map< unsigned int, edm::PSimHitContainer > superlayer_to_hits_
double p1[4]
Definition: TauolaWrapper.h:89
T get() const
Definition: EventSetup.h:71
float simHitsMeanWire(const edm::PSimHitContainer &sim_hits) const
std::vector< PSimHit > PSimHitContainer
void init(const edm::Event &e, const edm::EventSetup &eventSetup)
initialize the event
void init(const edm::Event &e, const edm::EventSetup &eventSetup)
initialize the event
int station() const
Return the station number.
Definition: DTChamberId.h:51
bool isValid() const
Definition: ESHandle.h:44
static int toDTType(int wh, int st)
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
const edm::PSimHitContainer & hitsInChamber(unsigned int) const
void dtChamberIdsToString(const std::set< unsigned int > &) const
def move(src, dest)
Definition: eostools.py:511
int nSuperLayersWithHitsInChamber(unsigned int) const
std::set< unsigned int > superlayerIds() const
int nCellsWithHitsInLayer(unsigned int) const
double p3[4]
Definition: TauolaWrapper.h:91
std::map< unsigned int, edm::PSimHitContainer > layer_to_hits_