CMS 3D CMS Logo

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