CMS 3D CMS Logo

GEMSegmentAlgorithm.cc
Go to the documentation of this file.
1 
8 #include "GEMSegmentAlgorithm.h"
9 #include "MuonSegFit.h"
10 
16 
17 #include <algorithm>
18 #include <cmath>
19 #include <iostream>
20 #include <string>
21 
22 /* Constructor
23  *
24  */
26  : GEMSegmentAlgorithmBase(ps), myName("GEMSegmentAlgorithm") {
27  minHitsPerSegment = ps.getParameter<unsigned int>("minHitsPerSegment");
28  preClustering = ps.getParameter<bool>("preClustering");
29  dXclusBoxMax = ps.getParameter<double>("dXclusBoxMax");
30  dYclusBoxMax = ps.getParameter<double>("dYclusBoxMax");
31  preClustering_useChaining = ps.getParameter<bool>("preClusteringUseChaining");
32  dPhiChainBoxMax = ps.getParameter<double>("dPhiChainBoxMax");
33  dEtaChainBoxMax = ps.getParameter<double>("dEtaChainBoxMax");
34  maxRecHitsInCluster = ps.getParameter<int>("maxRecHitsInCluster");
35  clusterOnlySameBXRecHits = ps.getParameter<bool>("clusterOnlySameBXRecHits");
36 
37  // maybe to be used in the future ???
38  // Pruning = ps.getParameter<bool>("Pruning");
39  // BrutePruning = ps.getParameter<bool>("BrutePruning");
40  // maxRecHitsInCluster is the maximal number of hits in a precluster that is being processed
41  // This cut is intended to remove messy events. Currently nothing is returned if there are
42  // more that maxRecHitsInCluster hits. It could be useful to return an estimate of the
43  // cluster position, which is available.
44  // maxRecHitsInCluster = ps.getParameter<int>("maxRecHitsInCluster");
45  // onlyBestSegment = ps.getParameter<bool>("onlyBestSegment");
46 
47  // CSC uses pruning to remove clearly bad hits, using as much info from the rechits as possible: charge, position, timing, ...
48  // In fits with bad chi^2 they look for the worst hit (hit with abnormally large residual)
49  // if worst hit was found, refit without worst hit and select if considerably better than original fit.
50 }
51 
52 /* Destructor
53  *
54  */
56 
57 std::vector<GEMSegment> GEMSegmentAlgorithm::run(const GEMEnsemble& ensemble, const EnsembleHitContainer& rechits) {
58  // pre-cluster rechits and loop over all sub clusters separately
59  std::vector<GEMSegment> segments_temp;
60  std::vector<GEMSegment> segments;
61  ProtoSegments rechits_clusters; // this is a collection of groups of rechits
62 
63  if (preClustering) {
64  // run a pre-clusterer on the given rechits to split obviously separated segment seeds:
66  // it uses X,Y,Z information; there are no configurable parameters used;
67  // the X, Y, Z "cuts" are just (much) wider than reasonable high pt segments
68  edm::LogVerbatim("GEMSegmentAlgorithm") << "[GEMSegmentAlgorithm::run] preClustering :: use Chaining";
69  rechits_clusters = this->chainHits(ensemble, rechits);
70  } else {
71  // it uses X,Y information + configurable parameters
72  edm::LogVerbatim("GEMSegmentAlgorithm") << "[GEMSegmentAlgorithm::run] Clustering";
73  rechits_clusters = this->clusterHits(ensemble, rechits);
74  }
75  // loop over the found clusters:
76  edm::LogVerbatim("GEMSegmentAlgorithm") << "[GEMSegmentAlgorithm::run] Loop over clusters and build segments";
77  for (auto sub_rechits = rechits_clusters.begin(); sub_rechits != rechits_clusters.end(); ++sub_rechits) {
78  // clear the buffer for the subset of segments:
79  segments_temp.clear();
80  // build the subset of segments:
81  this->buildSegments(ensemble, (*sub_rechits), segments_temp);
82  // add the found subset of segments to the collection of all segments in this chamber:
83  segments.insert(segments.end(), segments_temp.begin(), segments_temp.end());
84  }
85 
86  return segments;
87  } else {
88  this->buildSegments(ensemble, rechits, segments);
89  return segments;
90  }
91 }
92 
93 // ********************************************************************;
96  // think how to implement BX requirement here
97 
98  ProtoSegments rechits_clusters; // this is a collection of groups of rechits
99 
100  float dXclus_box = 0.0;
101  float dYclus_box = 0.0;
102 
104  seeds.reserve(rechits.size());
105 
106  std::vector<float> running_meanX;
107  running_meanX.reserve(rechits.size());
108  std::vector<float> running_meanY;
109  running_meanY.reserve(rechits.size());
110 
111  std::vector<float> seed_minX;
112  seed_minX.reserve(rechits.size());
113  std::vector<float> seed_maxX;
114  seed_maxX.reserve(rechits.size());
115  std::vector<float> seed_minY;
116  seed_minY.reserve(rechits.size());
117  std::vector<float> seed_maxY;
118  seed_maxY.reserve(rechits.size());
119 
120  // split rechits into subvectors and return vector of vectors:
121  // Loop over rechits
122  // Create one seed per hit
123  for (unsigned int i = 0; i < rechits.size(); ++i) {
124  seeds.push_back(EnsembleHitContainer(1, rechits[i]));
125 
126  GEMDetId rhID = rechits[i]->gemId();
127  const GEMEtaPartition* rhEP = (ensemble.second.find(rhID.rawId()))->second;
128  if (!rhEP)
129  throw cms::Exception("GEMEtaPartition not found")
130  << "Corresponding GEMEtaPartition to GEMDetId: " << rhID << " not found in the GEMEnsemble";
131  const GEMSuperChamber* rhCH = ensemble.first;
132  LocalPoint rhLP_inEtaPartFrame = rechits[i]->localPosition();
133  GlobalPoint rhGP_inCMSFrame = rhEP->toGlobal(rhLP_inEtaPartFrame);
134  LocalPoint rhLP_inChamberFrame = rhCH->toLocal(rhGP_inCMSFrame);
135 
136  running_meanX.push_back(rhLP_inChamberFrame.x());
137  running_meanY.push_back(rhLP_inChamberFrame.y());
138 
139  // set min/max X and Y for box containing the hits in the precluster:
140  seed_minX.push_back(rhLP_inChamberFrame.x());
141  seed_maxX.push_back(rhLP_inChamberFrame.x());
142  seed_minY.push_back(rhLP_inChamberFrame.y());
143  seed_maxY.push_back(rhLP_inChamberFrame.y());
144  }
145 
146  // merge clusters that are too close
147  // measure distance between final "running mean"
148  for (size_t NNN = 0; NNN < seeds.size(); ++NNN) {
149  for (size_t MMM = NNN + 1; MMM < seeds.size(); ++MMM) {
150  if (running_meanX[MMM] == running_max || running_meanX[NNN] == running_max) {
151  LogDebug("GEMSegmentAlgorithm") << "[GEMSegmentAlgorithm::clusterHits]: ALARM! Skipping used seeds, this "
152  "should not happen - inform developers!";
153  continue; //skip seeds that have been used
154  }
155 
156  // calculate cut criteria for simple running mean distance cut:
157  //dXclus = fabs(running_meanX[NNN] - running_meanX[MMM]);
158  //dYclus = fabs(running_meanY[NNN] - running_meanY[MMM]);
159  // calculate minmal distance between precluster boxes containing the hits:
160  if (running_meanX[NNN] > running_meanX[MMM])
161  dXclus_box = seed_minX[NNN] - seed_maxX[MMM];
162  else
163  dXclus_box = seed_minX[MMM] - seed_maxX[NNN];
164  if (running_meanY[NNN] > running_meanY[MMM])
165  dYclus_box = seed_minY[NNN] - seed_maxY[MMM];
166  else
167  dYclus_box = seed_minY[MMM] - seed_maxY[NNN];
168 
169  if (dXclus_box < dXclusBoxMax && dYclus_box < dYclusBoxMax) {
170  // merge clusters!
171  // merge by adding seed NNN to seed MMM and erasing seed NNN
172 
173  // calculate running mean for the merged seed:
174  if (seeds[NNN].size() + seeds[MMM].size() != 0) {
175  running_meanX[MMM] = (running_meanX[NNN] * seeds[NNN].size() + running_meanX[MMM] * seeds[MMM].size()) /
176  (seeds[NNN].size() + seeds[MMM].size());
177  running_meanY[MMM] = (running_meanY[NNN] * seeds[NNN].size() + running_meanY[MMM] * seeds[MMM].size()) /
178  (seeds[NNN].size() + seeds[MMM].size());
179  }
180 
181  // update min/max X and Y for box containing the hits in the merged cluster:
182  if (seed_minX[NNN] < seed_minX[MMM])
183  seed_minX[MMM] = seed_minX[NNN];
184  if (seed_maxX[NNN] > seed_maxX[MMM])
185  seed_maxX[MMM] = seed_maxX[NNN];
186  if (seed_minY[NNN] < seed_minY[MMM])
187  seed_minY[MMM] = seed_minY[NNN];
188  if (seed_maxY[NNN] > seed_maxY[MMM])
189  seed_maxY[MMM] = seed_maxY[NNN];
190 
191  // add seed NNN to MMM (lower to larger number)
192  seeds[MMM].insert(seeds[MMM].end(), seeds[NNN].begin(), seeds[NNN].end());
193 
194  // mark seed NNN as used (at the moment just set running mean to 999999.)
195  running_meanX[NNN] = running_max;
196  running_meanY[NNN] = running_max;
197  // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to
198  // next seed (NNN+1)
199  break;
200  }
201  }
202  }
203 
204  // hand over the final seeds to the output
205  // would be more elegant if we could do the above step with
206  // erasing the merged ones, rather than the
207  for (size_t NNN = 0; NNN < seeds.size(); ++NNN) {
208  if (running_meanX[NNN] == running_max)
209  continue; //skip seeds that have been marked as used up in merging
210  rechits_clusters.push_back(seeds[NNN]);
211  }
212 
213  return rechits_clusters;
214 }
215 
217  const EnsembleHitContainer& rechits) {
218  ProtoSegments rechits_chains;
220  seeds.reserve(rechits.size());
221  std::vector<bool> usedCluster(rechits.size(), false);
222 
223  // split rechits into subvectors and return vector of vectors:
224  // Loop over rechits
225  // Create one seed per hit
226  for (unsigned int i = 0; i < rechits.size(); ++i)
227  seeds.push_back(EnsembleHitContainer(1, rechits[i]));
228 
229  // merge chains that are too close ("touch" each other)
230  for (size_t NNN = 0; NNN < seeds.size(); ++NNN) {
231  for (size_t MMM = NNN + 1; MMM < seeds.size(); ++MMM) {
232  if (usedCluster[MMM] || usedCluster[NNN]) {
233  continue;
234  }
235  // all is in the way we define "good";
236  // try not to "cluster" the hits but to "chain" them;
237  // it does the clustering but also does a better job
238  // for inclined tracks (not clustering them together;
239  // crossed tracks would be still clustered together)
240  // 22.12.09: In fact it is not much more different
241  // than the "clustering", we just introduce another
242  // variable in the game - Z. And it makes sense
243  // to re-introduce Y (or actually wire group mumber)
244  // in a similar way as for the strip number - see
245  // the code below.
246  bool goodToMerge = isGoodToMerge(ensemble, seeds[NNN], seeds[MMM]);
247  if (goodToMerge) {
248  // merge chains!
249  // merge by adding seed NNN to seed MMM and erasing seed NNN
250 
251  // add seed NNN to MMM (lower to larger number)
252  seeds[MMM].insert(seeds[MMM].end(), seeds[NNN].begin(), seeds[NNN].end());
253 
254  // mark seed NNN as used
255  usedCluster[NNN] = true;
256  // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to
257  // next seed (NNN+1)
258  break;
259  }
260  }
261  }
262 
263  // hand over the final seeds to the output
264  // would be more elegant if we could do the above step with
265  // erasing the merged ones, rather than the
266  for (size_t NNN = 0; NNN < seeds.size(); ++NNN) {
267  if (usedCluster[NNN])
268  continue; //skip seeds that have been marked as used up in merging
269  rechits_chains.push_back(seeds[NNN]);
270  }
271 
272  //***************************************************************
273 
274  return rechits_chains;
275 }
276 
278  const EnsembleHitContainer& newChain,
279  const EnsembleHitContainer& oldChain) {
280  bool phiRequirementOK = false; // once it is true in the loop, it is ok to merge
281  bool etaRequirementOK = false; // once it is true in the loop, it is ok to merge
282  bool bxRequirementOK = false; // once it is true in the loop, it is ok to merge
283 
284  for (size_t iRH_new = 0; iRH_new < newChain.size(); ++iRH_new) {
285  int layer_new = (newChain[iRH_new]->gemId().station() - 1) * 2 + newChain[iRH_new]->gemId().layer();
286 
287  const GEMEtaPartition* rhEP = (ensemble.second.find(newChain[iRH_new]->gemId().rawId()))->second;
288  GlobalPoint pos_new = rhEP->toGlobal(newChain[iRH_new]->localPosition());
289 
290  for (size_t iRH_old = 0; iRH_old < oldChain.size(); ++iRH_old) {
291  int layer_old = (oldChain[iRH_old]->gemId().station() - 1) * 2 + oldChain[iRH_old]->gemId().layer();
292  // Layers - hits on the same layer should not be allowed ==> if abs(layer_new - layer_old) > 0 is ok. if = 0 is false
293  if (layer_new == layer_old)
294  return false;
295 
296  const GEMEtaPartition* oldrhEP = (ensemble.second.find(oldChain[iRH_old]->gemId().rawId()))->second;
297  GlobalPoint pos_old = oldrhEP->toGlobal(oldChain[iRH_old]->localPosition());
298 
299  // Eta & Phi- to be chained, two hits need also to be "close" in phi and eta
300  if (phiRequirementOK == false)
301  phiRequirementOK = std::abs(reco::deltaPhi(float(pos_new.phi()), float(pos_old.phi()))) < dPhiChainBoxMax;
302  if (etaRequirementOK == false)
303  etaRequirementOK = std::abs(pos_new.eta() - pos_old.eta()) < dEtaChainBoxMax;
304  // and they should have a time difference compatible with the hypothesis
305  // that the rechits originate from the same particle, but were detected in different layers
306  if (bxRequirementOK == false) {
308  bxRequirementOK = true;
309  } else {
310  if (newChain[iRH_new]->BunchX() == oldChain[iRH_old]->BunchX())
311  bxRequirementOK = true;
312  }
313  }
314 
315  if (phiRequirementOK && etaRequirementOK && bxRequirementOK)
316  return true;
317  }
318  }
319  return false;
320 }
321 
324  std::vector<GEMSegment>& gemsegs) {
325  if (rechits.size() < minHitsPerSegment)
326  return;
327 
329  proto_segment.clear();
330 
331  // select hits from the ensemble and sort it
332  const GEMSuperChamber* suCh = ensemble.first;
333  for (auto rh = rechits.begin(); rh != rechits.end(); rh++) {
334  proto_segment.push_back(*rh);
335 
336  // for segFit - using local point in chamber frame
337  const GEMEtaPartition* thePartition = (ensemble.second.find((*rh)->gemId()))->second;
338  GlobalPoint gp = thePartition->toGlobal((*rh)->localPosition());
339  const LocalPoint lp = suCh->toLocal(gp);
340 
341  GEMRecHit* newRH = (*rh)->clone();
342  newRH->setPosition(lp);
343  MuonSegFit::MuonRecHitPtr trkRecHit(newRH);
344  muonRecHits.push_back(trkRecHit);
345  }
346 
347 #ifdef EDM_ML_DEBUG // have lines below only compiled when in debug mode
348  edm::LogVerbatim("GEMSegmentAlgorithm")
349  << "[GEMSegmentAlgorithm::buildSegments] will now try to fit a GEMSegment from collection of " << rechits.size()
350  << " GEM RecHits";
351  for (auto rh = rechits.begin(); rh != rechits.end(); ++rh) {
352  auto gemid = (*rh)->gemId();
353  auto rhLP = (*rh)->localPosition();
354  edm::LogVerbatim("GEMSegmentAlgorithm")
355  << "[RecHit :: Loc x = " << std::showpos << std::setw(9) << rhLP.x() << " Loc y = " << std::showpos
356  << std::setw(9) << rhLP.y() << " BX = " << std::showpos << (*rh)->BunchX() << " -- " << gemid.rawId() << " = "
357  << gemid << " ]";
358  }
359 #endif
360 
361  // The actual fit on all hits of the vector of the selected Tracking RecHits:
362  sfit_ = std::make_unique<MuonSegFit>(muonRecHits);
363  bool goodfit = sfit_->fit();
364  edm::LogVerbatim("GEMSegmentAlgorithm")
365  << "[GEMSegmentAlgorithm::buildSegments] GEMSegment fit done :: fit is good = " << goodfit;
366 
367  // quit function if fit was not OK
368  if (!goodfit) {
369  for (auto rh : muonRecHits)
370  rh.reset();
371  return;
372  }
373  // obtain all information necessary to make the segment:
374  LocalPoint protoIntercept = sfit_->intercept();
375  LocalVector protoDirection = sfit_->localdir();
376  AlgebraicSymMatrix protoErrors = sfit_->covarianceMatrix();
377  double protoChi2 = sfit_->chi2();
378 
379  // Calculate the bunch crossing of the GEM Segment
380  float bx = 0.0;
381  for (auto rh = rechits.begin(); rh != rechits.end(); ++rh) {
382  bx += (*rh)->BunchX();
383  }
384  if (!rechits.empty())
385  bx = bx * 1.0 / (rechits.size());
386 
387  // Calculate the central value and uncertainty of the segment time
388  // if we want to study impact of 2-3ns time resolution on GEM Segment
389  // (if there will be TDCs in readout and not just BX determination)
390  // then implement tof() method for rechits and use this here
391  /*`
392  float averageTime=0.;
393  for (auto rh=rechits.begin(); rh!=rechits.end(); ++rh){
394  GEMEtaPartition thePartition = ensemble.second.find((*rh)->gemId));
395  GlobalPoint pos = (thePartition->toGlobal((*rh)->localPosition());
396  float tof = pos.mag() * 0.01 / 0.2997925 + 25.0*(*rh)->BunchX();
397  averageTime += pos;
398  }
399  if(rechits.size() != 0) averageTime=averageTime/(rechits.size());
400  float timeUncrt=0.;
401  for (auto rh=rechits.begin(); rh!=rechits.end(); ++rh){
402  GEMEtaPartition thePartition = ensemble.second.find((*rh)->gemId));
403  GlobalPoint pos = (thePartition->toGlobal((*rh)->localPosition());
404  float tof = pos.mag() * 0.01 / 0.2997925 + 25.0*(*rh)->BunchX();
405  timeUncrt += pow(tof-averageTime,2);
406  }
407  if(rechits.size() != 0) timeUncrt=timeUncrt/(rechits.size());
408  timeUncrt = sqrt(timeUncrt);
409  */
410 
411  // save all information inside GEMCSCSegment
412  edm::LogVerbatim("GEMSegmentAlgorithm")
413  << "[GEMSegmentAlgorithm::buildSegments] will now wrap fit info in GEMSegment dataformat";
414  GEMSegment tmp(proto_segment, protoIntercept, protoDirection, protoErrors, protoChi2, bx);
415  // GEMSegment tmp(proto_segment, protoIntercept, protoDirection, protoErrors, protoChi2, averageTime, timeUncrt);
416 
417  edm::LogVerbatim("GEMSegmentAlgorithm")
418  << "[GEMSegmentAlgorithm::buildSegments] GEMSegment made in " << tmp.gemDetId();
419  edm::LogVerbatim("GEMSegmentAlgorithm") << "[GEMSegmentAlgorithm::buildSegments] " << tmp;
420 
421  for (auto rh : muonRecHits)
422  rh.reset();
423  gemsegs.push_back(tmp);
424 }
size
Write out results.
std::shared_ptr< TrackingRecHit > MuonRecHitPtr
Definition: MuonSegFit.h:38
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
~GEMSegmentAlgorithm() override
Destructor.
ProtoSegments chainHits(const GEMEnsemble &ensemble, const EnsembleHitContainer &rechits)
T eta() const
Definition: PV3DBase.h:73
std::vector< const GEMRecHit * > EnsembleHitContainer
Typedefs.
std::pair< const GEMSuperChamber *, std::map< uint32_t, const GEMEtaPartition * > > GEMEnsemble
U second(std::pair< T, U > const &p)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
bool isGoodToMerge(const GEMEnsemble &ensemble, const EnsembleHitContainer &newChain, const EnsembleHitContainer &oldChain)
unsigned int minHitsPerSegment
GEMRecHit * clone() const override
Definition: GEMRecHit.cc:58
std::vector< MuonRecHitPtr > MuonRecHitContainer
Definition: MuonSegFit.h:39
ProtoSegments clusterHits(const GEMEnsemble &ensemble, const EnsembleHitContainer &rechits)
Utility functions.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< EnsembleHitContainer > ProtoSegments
std::unique_ptr< MuonSegFit > sfit_
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
GEMSegmentAlgorithm(const edm::ParameterSet &ps)
Constructor.
static constexpr float running_max
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
CLHEP::HepSymMatrix AlgebraicSymMatrix
EnsembleHitContainer proto_segment
void setPosition(LocalPoint pos)
Set local position.
Definition: GEMRecHit.h:53
tmp
align.sh
Definition: createJobs.py:716
std::vector< GEMSegment > run(const GEMEnsemble &ensemble, const EnsembleHitContainer &rechits) override
void buildSegments(const GEMEnsemble &ensemble, const EnsembleHitContainer &rechits, std::vector< GEMSegment > &gemsegs)
#define LogDebug(id)