CMS 3D CMS Logo

ME0SegmentAlgorithm.cc
Go to the documentation of this file.
1 
8 #include "ME0SegmentAlgorithm.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 {
27  debug = ps.getUntrackedParameter<bool>("ME0Debug");
28  minHitsPerSegment = ps.getParameter<unsigned int>("minHitsPerSegment");
29  preClustering = ps.getParameter<bool>("preClustering");
30  dXclusBoxMax = ps.getParameter<double>("dXclusBoxMax");
31  dYclusBoxMax = ps.getParameter<double>("dYclusBoxMax");
32  preClustering_useChaining = ps.getParameter<bool>("preClusteringUseChaining");
33  dPhiChainBoxMax = ps.getParameter<double>("dPhiChainBoxMax");
34  dEtaChainBoxMax = ps.getParameter<double>("dEtaChainBoxMax");
35  dTimeChainBoxMax = ps.getParameter<double>("dTimeChainBoxMax");
36  maxRecHitsInCluster = ps.getParameter<int>("maxRecHitsInCluster");
37 
38  edm::LogVerbatim("ME0SegmentAlgorithm") << "[ME0SegmentAlgorithm::ctor] Parameters to build segments :: "
39  << "preClustering = "<<preClustering<<" preClustering_useChaining = "<<preClustering_useChaining
40  <<" dPhiChainBoxMax = "<<dPhiChainBoxMax<<" dEtaChainBoxMax = "<<dEtaChainBoxMax<<" dTimeChainBoxMax = "<<dTimeChainBoxMax
41  <<" minHitsPerSegment = "<<minHitsPerSegment<<" maxRecHitsInCluster = "<<maxRecHitsInCluster;
42 }
43 
44 /* Destructor
45  *
46  */
48 }
49 
50 
52 
53  #ifdef EDM_ML_DEBUG // have lines below only compiled when in debug mode
54  ME0DetId chId(chamber->id());
55  edm::LogVerbatim("ME0SegAlgoMM") << "[ME0SegmentAlgorithm::run] build segments in chamber " << chId << " which contains "<<rechits.size()<<" rechits";
56  for (auto rh=rechits.begin(); rh!=rechits.end(); ++rh){
57  auto me0id = rh->rh->me0Id();
58  auto rhLP = rh->lp;
59  edm::LogVerbatim("ME0SegmentAlgorithm") << "[RecHit :: Loc x = "<<std::showpos<<std::setw(9)<<rhLP.x()<<" Loc y = "<<std::showpos<<std::setw(9)<<rhLP.y()
60  <<" Time = "<<std::showpos<<rh->rh->tof()<<" -- "<<me0id.rawId()<<" = "<<me0id<<" ]";
61  }
62  #endif
63 
64  // pre-cluster rechits and loop over all sub clusters separately
65  std::vector<ME0Segment> segments_temp;
66  std::vector<ME0Segment> segments;
67  ProtoSegments rechits_clusters; // this is a collection of groups of rechits
68 
69  if(preClustering) {
70  // run a pre-clusterer on the given rechits to split obviously separated segment seeds:
72  // it uses X,Y,Z information; there are no configurable parameters used;
73  // the X, Y, Z "cuts" are just (much) wider than reasonable high pt segments
74  rechits_clusters = this->chainHits(chamber, rechits );
75  }
76  else{
77  // it uses X,Y information + configurable parameters
78  rechits_clusters = this->clusterHits(rechits );
79  }
80  // loop over the found clusters:
81  for(auto sub_rechits = rechits_clusters.begin(); sub_rechits != rechits_clusters.end(); ++sub_rechits ) {
82  // clear the buffer for the subset of segments:
83  segments_temp.clear();
84  // build the subset of segments:
85  this->buildSegments(chamber,(*sub_rechits), segments_temp );
86  // add the found subset of segments to the collection of all segments in this chamber:
87  segments.insert( segments.end(), segments_temp.begin(), segments_temp.end() );
88  }
89  return segments;
90  }
91  else {
92  HitAndPositionPtrContainer ptrRH; ptrRH.reserve(rechits.size());
93  for(const auto& rh : rechits) ptrRH.push_back(&rh);
94  this->buildSegments(chamber,ptrRH, segments);
95  return segments;
96  }
97 }
98 
99 
100 // ********************************************************************;
103 
104  ProtoSegments rechits_clusters; // this is a collection of groups of rechits
105 
106  float dXclus_box = 0.0;
107  float dYclus_box = 0.0;
108 
109  ProtoSegments seeds; seeds.reserve(rechits.size());
110 
111  std::vector<float> running_meanX; running_meanX.reserve(rechits.size());
112  std::vector<float> running_meanY; running_meanY.reserve(rechits.size());
113 
114  std::vector<float> seed_minX; seed_minX.reserve(rechits.size());
115  std::vector<float> seed_maxX; seed_maxX.reserve(rechits.size());
116  std::vector<float> seed_minY; seed_minY.reserve(rechits.size());
117  std::vector<float> seed_maxY; seed_maxY.reserve(rechits.size());
118 
119  // split rechits into subvectors and return vector of vectors:
120  // Loop over rechits
121  // Create one seed per hit
122  for(unsigned int i = 0; i < rechits.size(); ++i) {
123  seeds.push_back(HitAndPositionPtrContainer(1,&rechits[i]));
124 
125  // First added hit in seed defines the mean to which the next hit is compared
126  // for this seed.
127 
128  running_meanX.push_back( rechits[i].lp.x() );
129  running_meanY.push_back( rechits[i].lp.y() );
130 
131  // set min/max X and Y for box containing the hits in the precluster:
132  seed_minX.push_back( rechits[i].lp.x() );
133  seed_maxX.push_back( rechits[i].lp.x() );
134  seed_minY.push_back( rechits[i].lp.y() );
135  seed_maxY.push_back( rechits[i].lp.y() );
136  }
137 
138  // merge clusters that are too close
139  // measure distance between final "running mean"
140  for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
141  for(size_t MMM = NNN+1; MMM < seeds.size(); ++MMM) {
142  if(running_meanX[MMM] == running_max || running_meanX[NNN] == running_max ) {
143  LogDebug("ME0SegmentAlgorithm") << "[ME0SegmentAlgorithm::clusterHits]: ALARM! Skipping used seeds, this should not happen - inform developers!";
144  continue; //skip seeds that have been used
145  }
146 
147  // calculate cut criteria for simple running mean distance cut:
148  //dXclus = fabs(running_meanX[NNN] - running_meanX[MMM]);
149  //dYclus = fabs(running_meanY[NNN] - running_meanY[MMM]);
150  // calculate minmal distance between precluster boxes containing the hits:
151  if ( running_meanX[NNN] > running_meanX[MMM] ) dXclus_box = seed_minX[NNN] - seed_maxX[MMM];
152  else dXclus_box = seed_minX[MMM] - seed_maxX[NNN];
153  if ( running_meanY[NNN] > running_meanY[MMM] ) dYclus_box = seed_minY[NNN] - seed_maxY[MMM];
154  else dYclus_box = seed_minY[MMM] - seed_maxY[NNN];
155 
156 
157  if( dXclus_box < dXclusBoxMax && dYclus_box < dYclusBoxMax ) {
158  // merge clusters!
159  // merge by adding seed NNN to seed MMM and erasing seed NNN
160 
161  // calculate running mean for the merged seed:
162  if(seeds[NNN].size()+seeds[MMM].size() != 0) {
163  running_meanX[MMM] = (running_meanX[NNN]*seeds[NNN].size() + running_meanX[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
164  running_meanY[MMM] = (running_meanY[NNN]*seeds[NNN].size() + running_meanY[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
165  }
166 
167  // update min/max X and Y for box containing the hits in the merged cluster:
168  if ( seed_minX[NNN] < seed_minX[MMM] ) seed_minX[MMM] = seed_minX[NNN];
169  if ( seed_maxX[NNN] > seed_maxX[MMM] ) seed_maxX[MMM] = seed_maxX[NNN];
170  if ( seed_minY[NNN] < seed_minY[MMM] ) seed_minY[MMM] = seed_minY[NNN];
171  if ( seed_maxY[NNN] > seed_maxY[MMM] ) seed_maxY[MMM] = seed_maxY[NNN];
172 
173  // add seed NNN to MMM (lower to larger number)
174  seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
175 
176  // mark seed NNN as used (at the moment just set running mean to 999999.)
177  running_meanX[NNN] = running_max;
178  running_meanY[NNN] = running_max;
179  // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to
180  // next seed (NNN+1)
181  break;
182  }
183  }
184  }
185 
186  // hand over the final seeds to the output
187  // would be more elegant if we could do the above step with
188  // erasing the merged ones, rather than the
189  for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
190  if(running_meanX[NNN] == running_max) continue; //skip seeds that have been marked as used up in merging
191  rechits_clusters.push_back(seeds[NNN]);
192  }
193 
194  return rechits_clusters;
195 }
196 
197 
200 
201  ProtoSegments rechits_chains;
202  ProtoSegments seeds;
203  seeds.reserve(rechits.size());
204  std::vector<bool> usedCluster(rechits.size(),false);
205 
206  // split rechits into subvectors and return vector of vectors:
207  // Loop over rechits
208  // Create one seed per hit
209  for ( unsigned int i=0; i<rechits.size(); ++i){
210  if(std::abs(rechits[i].rh->tof()) > dTimeChainBoxMax) continue;
211  seeds.push_back(HitAndPositionPtrContainer(1,&rechits[i]));
212  }
213 
214  // merge chains that are too close ("touch" each other)
215  for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
216  for(size_t MMM = NNN+1; MMM < seeds.size(); ++MMM) {
217  if(usedCluster[MMM] || usedCluster[NNN]){
218  continue;
219  }
220  // all is in the way we define "good";
221  // try not to "cluster" the hits but to "chain" them;
222  // it does the clustering but also does a better job
223  // for inclined tracks (not clustering them together;
224  // crossed tracks would be still clustered together)
225  // 22.12.09: In fact it is not much more different
226  // than the "clustering", we just introduce another
227  // variable in the game - Z. And it makes sense
228  // to re-introduce Y (or actually wire group mumber)
229  // in a similar way as for the strip number - see
230  // the code below.
231  bool goodToMerge = isGoodToMerge(chamber, seeds[NNN], seeds[MMM]);
232  if(goodToMerge){
233  // merge chains!
234  // merge by adding seed NNN to seed MMM and erasing seed NNN
235 
236  // add seed NNN to MMM (lower to larger number)
237  seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
238 
239  // mark seed NNN as used
240  usedCluster[NNN] = true;
241  // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to
242  // next seed (NNN+1)
243  break;
244  }
245 
246  }
247  }
248 
249  // hand over the final seeds to the output
250  // would be more elegant if we could do the above step with
251  // erasing the merged ones, rather than the
252 
253  for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
254  if(usedCluster[NNN]) continue; //skip seeds that have been marked as used up in merging
255  rechits_chains.push_back(seeds[NNN]);
256  }
257 
258  //***************************************************************
259 
260  return rechits_chains;
261 }
262 
264  for(size_t iRH_new = 0;iRH_new<newChain.size();++iRH_new){
265  GlobalPoint pos_new = newChain[iRH_new]->gp;
266 
267  for(size_t iRH_old = 0;iRH_old<oldChain.size();++iRH_old){
268  GlobalPoint pos_old = oldChain[iRH_old]->gp;
269  // to be chained, two hits need to be in neighbouring layers...
270  // or better allow few missing layers (upto 3 to avoid inefficiencies);
271  // however we'll not make an angle correction because it
272  // worsen the situation in some of the "regular" cases
273  // (not making the correction means that the conditions for
274  // forming a cluster are different if we have missing layers -
275  // this could affect events at the boundaries )
276 
277  // to be chained, two hits need also to be "close" in phi and eta
278  if (std::abs(reco::deltaPhi( float(pos_new.phi()), float(pos_old.phi()) )) >= dPhiChainBoxMax) continue;
279  if (std::abs(pos_new.eta()-pos_old.eta()) >= dEtaChainBoxMax) continue;
280  // and the difference in layer index should be < (nlayers-1)
281  if (std::abs(int(newChain[iRH_new]->layer) - int(oldChain[iRH_old]->layer)) >= (chamber->id().nlayers()-1)) continue;
282  // and they should have a time difference compatible with the hypothesis
283  // that the rechits originate from the same particle, but were detected in different layers
284  if (std::abs(newChain[iRH_new]->rh->tof() - oldChain[iRH_old]->rh->tof()) >= dTimeChainBoxMax) continue;
285 
286  return true;
287  }
288  }
289  return false;
290 }
291 
292 void ME0SegmentAlgorithm::buildSegments(const ME0Chamber * chamber, const HitAndPositionPtrContainer& rechits, std::vector<ME0Segment>& me0segs) {
293  if (rechits.size() < minHitsPerSegment) return;
294 
295 #ifdef EDM_ML_DEBUG // have lines below only compiled when in debug mode
296  edm::LogVerbatim("ME0SegmentAlgorithm") << "[ME0SegmentAlgorithm::buildSegments] will now try to fit a ME0Segment from collection of "<<rechits.size()<<" ME0 RecHits";
297  for (auto rh=rechits.begin(); rh!=rechits.end(); ++rh){
298  auto me0id = (*rh)->rh->me0Id();
299  auto rhLP = (*rh)->lp;
300  edm::LogVerbatim("ME0SegmentAlgorithm") << "[RecHit :: Loc x = "<<std::showpos<<std::setw(9)<<rhLP.x()<<" Loc y = "<<std::showpos<<std::setw(9)<<rhLP.y()<<" Time = "<<std::showpos<<(*rh)->rh->tof()<<" -- "<<me0id.rawId()<<" = "<<me0id<<" ]";
301  }
302 #endif
303 
305  std::vector<const ME0RecHit*> bareRHs;
306 
307  // select hits from the ensemble and sort it
308  for (auto rh=rechits.begin(); rh!=rechits.end();rh++){
309  bareRHs.push_back((*rh)->rh );
310  // for segFit - using local point in chamber frame
311  ME0RecHit *newRH = (*rh)->rh->clone();
312  newRH->setPosition((*rh)->lp);
313 
314  MuonSegFit::MuonRecHitPtr trkRecHit(newRH);
315  muonRecHits.push_back(trkRecHit);
316  }
317 
318  // The actual fit on all hits of the vector of the selected Tracking RecHits:
319  sfit_ = std::make_unique<MuonSegFit>(muonRecHits);
320  bool goodfit = sfit_->fit();
321  edm::LogVerbatim("ME0SegmentAlgorithm") << "[ME0SegmentAlgorithm::buildSegments] ME0Segment fit done";
322 
323  // quit function if fit was not OK
324  if(!goodfit){
325  for (auto rh:muonRecHits) rh.reset();
326  return;
327  }
328 
329  // obtain all information necessary to make the segment:
330  LocalPoint protoIntercept = sfit_->intercept();
331  LocalVector protoDirection = sfit_->localdir();
332  AlgebraicSymMatrix protoErrors = sfit_->covarianceMatrix();
333  double protoChi2 = sfit_->chi2();
334  // Calculate the central value and uncertainty of the segment time
335  float averageTime=0.;
336  for (auto rh=rechits.begin(); rh!=rechits.end(); ++rh){
337  averageTime += (*rh)->rh->tof();
338  }
339  if(!rechits.empty()) averageTime=averageTime/(rechits.size());
340  float timeUncrt=0.;
341  for (auto rh=rechits.begin(); rh!=rechits.end(); ++rh){
342  timeUncrt += pow((*rh)->rh->tof()-averageTime,2);
343  }
344 
345  if(rechits.size() > 1) timeUncrt=timeUncrt/(rechits.size()-1);
346  timeUncrt = sqrt(timeUncrt);
347 
348  const float dPhi = chamber->computeDeltaPhi(protoIntercept,protoDirection);
349 
350  // save all information inside GEMCSCSegment
351  edm::LogVerbatim("ME0SegmentAlgorithm") << "[ME0SegmentAlgorithm::buildSegments] will now try to make ME0Segment from collection of "<<rechits.size()<<" ME0 RecHits";
352  ME0Segment tmp(bareRHs, protoIntercept, protoDirection, protoErrors, protoChi2, averageTime, timeUncrt, dPhi);
353 
354  edm::LogVerbatim("ME0SegmentAlgorithm") << "[ME0SegmentAlgorithm::buildSegments] ME0Segment made";
355  edm::LogVerbatim("ME0SegmentAlgorithm") << "[ME0SegmentAlgorithm::buildSegments] "<<tmp;
356 
357  for (auto rh:muonRecHits) rh.reset();
358  me0segs.push_back(tmp);
359  return;
360 }
#define LogDebug(id)
size
Write out results.
std::shared_ptr< TrackingRecHit > MuonRecHitPtr
Definition: MuonSegFit.h:41
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool isGoodToMerge(const ME0Chamber *chamber, const HitAndPositionPtrContainer &newChain, const HitAndPositionPtrContainer &oldChain)
void setPosition(LocalPoint pos)
Set local position.
Definition: ME0RecHit.h:72
void buildSegments(const ME0Chamber *chamber, const HitAndPositionPtrContainer &rechits, std::vector< ME0Segment > &me0segs)
unsigned int minHitsPerSegment
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
ME0DetId id() const
Return the ME0DetId of this chamber.
Definition: ME0Chamber.cc:14
std::vector< ME0Segment > run(const ME0Chamber *chamber, const HitAndPositionContainer &rechits) override
std::vector< HitAndPosition > HitAndPositionContainer
ProtoSegments chainHits(const ME0Chamber *chamber, const HitAndPositionContainer &rechits)
~ME0SegmentAlgorithm() override
Destructor.
std::vector< const HitAndPosition * > HitAndPositionPtrContainer
ProtoSegments clusterHits(const HitAndPositionContainer &rechits)
Utility functions.
float computeDeltaPhi(const LocalPoint &position, const LocalVector &direction) const
Definition: ME0Chamber.cc:82
std::unique_ptr< MuonSegFit > sfit_
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< MuonRecHitPtr > MuonRecHitContainer
Definition: MuonSegFit.h:42
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define end
Definition: vmac.h:39
std::vector< HitAndPositionPtrContainer > ProtoSegments
Typedefs.
ME0RecHit * clone() const override
Definition: ME0RecHit.cc:47
T eta() const
Definition: PV3DBase.h:76
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
int nlayers() const
For future modifications (implement more layers)
Definition: ME0DetId.h:82
#define begin
Definition: vmac.h:32
CLHEP::HepSymMatrix AlgebraicSymMatrix
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
ME0SegmentAlgorithm(const edm::ParameterSet &ps)
Constructor.