CMS 3D CMS Logo

DTCombinatorialPatternReco.cc
Go to the documentation of this file.
1 
7 /* This Class Header */
9 
10 /* Collaborating Class Header */
25 
26 /* C++ Headers */
27 #include <iterator>
28 using namespace std;
30 
31 /* ====================================================================== */
32 
35  : DTRecSegment2DBaseAlgo(pset), theAlgoName("DTCombinatorialPatternReco"), theDTGeometryToken(iCC.esConsumes()) {
36  theMaxAllowedHits = pset.getParameter<unsigned int>("MaxAllowedHits"); // 100
37  theAlphaMaxTheta = pset.getParameter<double>("AlphaMaxTheta"); // 0.1 ;
38  theAlphaMaxPhi = pset.getParameter<double>("AlphaMaxPhi"); // 1.0 ;
39  debug = pset.getUntrackedParameter<bool>("debug"); //true;
40  theUpdator = new DTSegmentUpdator(pset, iCC);
42  string theHitAlgoName = pset.getParameter<string>("recAlgo");
43  usePairs = !(theHitAlgoName == "DTNoDriftAlgo");
44 }
45 
48  delete theUpdator;
49  delete theCleaner;
50 }
51 
52 /* Operations */
54  const std::vector<DTRecHit1DPair>& pairs) {
56  vector<std::shared_ptr<DTHitPairForFit>> hitsForFit = initHits(sl, pairs);
57 
58  vector<DTSegmentCand*> candidates = buildSegments(sl, hitsForFit);
59 
60  vector<DTSegmentCand*>::const_iterator cand = candidates.begin();
61  while (cand < candidates.end()) {
62  DTSLRecSegment2D* segment = (**cand);
63 
64  theUpdator->update(segment, false);
65 
66  result.push_back(segment);
67 
68  if (debug) {
69  cout << "Reconstructed 2D segments " << result.back() << endl;
70  }
71 
72  delete *(cand++); // delete the candidate!
73  }
74 
75  return result;
76 }
77 
79  // Get the DT Geometry
82 }
83 
84 vector<std::shared_ptr<DTHitPairForFit>> DTCombinatorialPatternReco::initHits(const DTSuperLayer* sl,
85  const std::vector<DTRecHit1DPair>& hits) {
86  vector<std::shared_ptr<DTHitPairForFit>> result;
87  for (vector<DTRecHit1DPair>::const_iterator hit = hits.begin(); hit != hits.end(); ++hit) {
88  result.push_back(std::make_shared<DTHitPairForFit>(*hit, *sl, theDTGeometry));
89  }
90  return result;
91 }
92 
94  const DTSuperLayer* sl, const std::vector<std::shared_ptr<DTHitPairForFit>>& hits) {
95  // clear the patterns tried
96  if (debug) {
97  cout << "theTriedPattern.size is " << theTriedPattern.size() << "\n";
98  }
99  theTriedPattern.clear();
100 
101  typedef vector<std::shared_ptr<DTHitPairForFit>> hitCont;
102  typedef hitCont::const_iterator hitIter;
103  vector<DTSegmentCand*> result;
104 
105  if (debug) {
106  cout << "buildSegments: " << sl->id() << " nHits " << hits.size() << endl;
107  for (vector<std::shared_ptr<DTHitPairForFit>>::const_iterator hit = hits.begin(); hit != hits.end(); ++hit)
108  cout << **hit << endl;
109  }
110 
111  // 10-Mar-2004 SL
112  // put a protection against heavily populated chambers, for which the segment
113  // building could lead to infinite memory usage...
114  if (hits.size() > theMaxAllowedHits) {
115  if (debug) {
116  cout << "Warning: this SuperLayer " << sl->id() << " has too many hits : " << hits.size() << " max allowed is "
117  << theMaxAllowedHits << endl;
118  cout << "Skipping segment reconstruction... " << endl;
119  }
120  return result;
121  }
122 
124  // compatible with them
125  for (hitCont::const_iterator firstHit = hits.begin(); firstHit != hits.end(); ++firstHit) {
126  for (hitCont::const_reverse_iterator lastHit = hits.rbegin(); (*lastHit) != (*firstHit); ++lastHit) {
127  //if ( (*lastHit)->id().layerId() == (*firstHit)->id().layerId() ) continue; // hits must be in different layers!
128  // hits must nor in the same nor in adiacent layers
129  if (fabs((*lastHit)->id().layerId() - (*firstHit)->id().layerId()) <= 1)
130  continue;
131  if (debug) {
132  cout << "Selected these two hits pair " << endl;
133  cout << "First " << *(*firstHit) << " Layer Id: " << (*firstHit)->id().layerId() << endl;
134  cout << "Last " << *(*lastHit) << " Layer Id: " << (*lastHit)->id().layerId() << endl;
135  }
136 
137  GlobalPoint IP;
138  float DAlphaMax;
139  if ((sl->id()).superlayer() == 2) // Theta SL
140  DAlphaMax = theAlphaMaxTheta;
141  else // Phi SL
142  DAlphaMax = theAlphaMaxPhi;
143 
145  for (int firstLR = 0; firstLR < 2; ++firstLR) {
146  for (int lastLR = 0; lastLR < 2; ++lastLR) {
147  // TODO move the global transformation in the DTHitPairForFit class
148  // when it will be moved I will able to remove the sl from the input parameter
149  GlobalPoint gposFirst = sl->toGlobal((*firstHit)->localPosition(codes[firstLR]));
150  GlobalPoint gposLast = sl->toGlobal((*lastHit)->localPosition(codes[lastLR]));
151 
152  GlobalVector gvec = gposLast - gposFirst;
153  GlobalVector gvecIP = gposLast - IP;
154 
155  // difference in angle measured
156  float DAlpha = fabs(gvec.theta() - gvecIP.theta());
157 
158  // cout << "DAlpha " << DAlpha << endl;
159  if (DAlpha < DAlphaMax) {
160  // create a segment hypotesis
161  // I don't need a true segment, just direction and position
162  LocalPoint posIni = (*firstHit)->localPosition(codes[firstLR]);
163  LocalVector dirIni = ((*lastHit)->localPosition(codes[lastLR]) - posIni).unit();
164 
165  // search for other compatible hits, with or without the L/R solved
166  vector<DTSegmentCand::AssPoint> assHits = findCompatibleHits(posIni, dirIni, hits);
167  if (debug)
168  cout << "compatible hits " << assHits.size() << endl;
169 
170  // get the best segment with these hits: it's just one!
171  // (is it correct?)
172  DTSegmentCand* seg = buildBestSegment(assHits, sl);
173 
174  if (seg) {
175  if (debug)
176  cout << "segment " << *seg << endl;
177 
178  // check if the chi2 and #hits are ok
179  if (!seg->good()) {
180  delete seg;
181  } else {
182  // remove duplicated segments (I know, would be better to do it before the
183  // fit...)
184  if (checkDoubleCandidates(result, seg)) {
185  // add to the vector of hypotesis
186  result.push_back(seg);
187  if (debug)
188  cout << "result is now " << result.size() << endl;
189  } else { // delete it!
190  delete seg;
191  if (debug)
192  cout << "already existing" << endl;
193  }
194  }
195  }
196  }
197  }
198  }
199  }
200  }
201  if (debug) {
202  for (vector<DTSegmentCand*>::const_iterator seg = result.begin(); seg != result.end(); ++seg)
203  cout << *(*seg) << endl;
204  }
205 
206  // now I have a couple of segment hypotesis, should check for ghost
208  if (debug) {
209  cout << "result no ghost " << result.size() << endl;
210  for (vector<DTSegmentCand*>::const_iterator seg = result.begin(); seg != result.end(); ++seg)
211  cout << *(*seg) << endl;
212  }
213 
214  return result;
215 }
216 
217 vector<DTSegmentCand::AssPoint> DTCombinatorialPatternReco::findCompatibleHits(
218  const LocalPoint& posIni, const LocalVector& dirIni, const vector<std::shared_ptr<DTHitPairForFit>>& hits) {
219  if (debug)
220  cout << "Pos: " << posIni << " Dir: " << dirIni << endl;
221  vector<DTSegmentCand::AssPoint> result;
222 
223  // counter to early-avoid double counting in hits pattern
225  int nCompatibleHits = 0;
226 
227  typedef vector<std::shared_ptr<DTHitPairForFit>> hitCont;
228  typedef hitCont::const_iterator hitIter;
229  for (hitIter hit = hits.begin(); hit != hits.end(); ++hit) {
230  pair<bool, bool> isCompatible = (*hit)->isCompatible(posIni, dirIni);
231  if (debug)
232  cout << "isCompatible " << isCompatible.first << " " << isCompatible.second << endl;
233 
234  // if only one of the two is compatible, then the LR is assigned,
235  // otherwise is undefined
236 
237  DTEnums::DTCellSide lrcode;
238  if (isCompatible.first && isCompatible.second) {
239  usePairs ? lrcode = DTEnums::undefLR : lrcode = DTEnums::Left; // if not usePairs then only use single side
240  tried.push_back(3);
241  nCompatibleHits++;
242  } else if (isCompatible.first) {
243  lrcode = DTEnums::Left;
244  tried.push_back(2);
245  nCompatibleHits++;
246  } else if (isCompatible.second) {
247  lrcode = DTEnums::Right;
248  tried.push_back(1);
249  nCompatibleHits++;
250  } else {
251  tried.push_back(0);
252  continue; // neither is compatible
253  }
254  result.push_back(DTSegmentCand::AssPoint(*hit, lrcode));
255  }
256 
257  // check if too few associated hits or pattern already tried
258  // TODO: two different if for nCompatibleHits < 3 =>printout and find ->
259  // printour
260  if (nCompatibleHits < 3) {
261  if (debug) {
262  cout << "Less than 3 hits ";
263  copy(tried.begin(), tried.end(), ostream_iterator<int>(std::cout));
264  cout << endl;
265  }
266  // No need to insert into the list of patterns, as you will never search for it.
267  //theTriedPattern.insert(tried);
268  } else {
269  // try to insert, return bool if already there
270  bool isnew = theTriedPattern.insert(tried).second;
271  if (isnew) {
272  if (debug) {
273  cout << "NOT Already tried ";
274  copy(tried.begin(), tried.end(), ostream_iterator<int>(std::cout));
275  cout << endl;
276  }
277  } else {
278  if (debug) {
279  cout << "Already tried ";
280  copy(tried.begin(), tried.end(), ostream_iterator<int>(std::cout));
281  cout << endl;
282  }
283  // empty the result vector
284  result.clear();
285  }
286  }
287  return result;
288 }
289 
290 DTSegmentCand* DTCombinatorialPatternReco::buildBestSegment(std::vector<DTSegmentCand::AssPoint>& hits,
291  const DTSuperLayer* sl) {
292  if (hits.size() < 3) {
293  //cout << "buildBestSegment: hits " << hits.size()<< endl;
294  return nullptr; // a least 3 point
295  }
296 
297  // hits with defined LR
298  vector<DTSegmentCand::AssPoint> points;
299 
300  // without: I store both L and R, a deque since I need front insertion and
301  // deletion
302  deque<std::shared_ptr<DTHitPairForFit>> pointsNoLR;
303 
304  // first add only the hits with LR assigned
305  for (vector<DTSegmentCand::AssPoint>::const_iterator hit = hits.begin(); hit != hits.end(); ++hit) {
306  if ((*hit).second != DTEnums::undefLR) {
307  points.push_back(*hit);
308  } else { // then also for the undef'd one
309  pointsNoLR.push_back((*hit).first);
310  }
311  }
312 
313  if (debug) {
314  cout << "points " << points.size() << endl;
315  cout << "pointsNoLR " << pointsNoLR.size() << endl;
316  }
317 
318  // build all possible candidates using L/R ambiguity
319  vector<DTSegmentCand*> candidates;
320 
321  buildPointsCollection(points, pointsNoLR, candidates, sl);
322 
323  if (debug)
324  cout << "candidates " << candidates.size() << endl;
325 
326  // so now I have build a given number of segments, I should find the best one,
327  // by #hits and chi2.
328  vector<DTSegmentCand*>::const_iterator bestCandIter = candidates.end();
329  double minChi2 = 999999.;
330  unsigned int maxNumHits = 0;
331  for (vector<DTSegmentCand*>::const_iterator iter = candidates.begin(); iter != candidates.end(); ++iter) {
332  if ((*iter)->nHits() == maxNumHits && (*iter)->chi2() < minChi2) {
333  minChi2 = (*iter)->chi2();
334  bestCandIter = iter;
335  } else if ((*iter)->nHits() > maxNumHits) {
336  maxNumHits = (*iter)->nHits();
337  minChi2 = (*iter)->chi2();
338  bestCandIter = iter;
339  }
340  }
341 
342  // delete all candidates but the best one!
343  for (vector<DTSegmentCand*>::iterator iter = candidates.begin(); iter != candidates.end(); ++iter)
344  if (iter != bestCandIter)
345  delete *iter;
346 
347  // return the best candate if any
348  if (bestCandIter != candidates.end()) {
349  return (*bestCandIter);
350  }
351  return nullptr;
352 }
353 
354 void DTCombinatorialPatternReco::buildPointsCollection(vector<DTSegmentCand::AssPoint>& points,
355  deque<std::shared_ptr<DTHitPairForFit>>& pointsNoLR,
356  vector<DTSegmentCand*>& candidates,
357  const DTSuperLayer* sl) {
358  if (debug) {
359  cout << "buildPointsCollection " << endl;
360  cout << "points: " << points.size() << " NOLR: " << pointsNoLR.size() << endl;
361  }
362  if (!pointsNoLR.empty()) { // still unassociated points!
363  std::shared_ptr<DTHitPairForFit> unassHit = pointsNoLR.front();
364  // try with the right
365  if (debug)
366  cout << "Right hit" << endl;
367  points.push_back(DTSegmentCand::AssPoint(unassHit, DTEnums::Right));
368  pointsNoLR.pop_front();
369  buildPointsCollection(points, pointsNoLR, candidates, sl);
370  pointsNoLR.push_front((unassHit));
371  points.pop_back();
372 
373  // try with the left
374  if (debug)
375  cout << "Left hit" << endl;
376  points.push_back(DTSegmentCand::AssPoint(unassHit, DTEnums::Left));
377  pointsNoLR.pop_front();
378  buildPointsCollection(points, pointsNoLR, candidates, sl);
379  pointsNoLR.push_front((unassHit));
380  points.pop_back();
381  } else { // all associated
382 
383  if (debug) {
384  cout << "The Hits were" << endl;
385  copy(points.begin(), points.end(), ostream_iterator<DTSegmentCand::AssPoint>(std::cout));
386  cout << "----" << endl;
387  cout << "All associated " << endl;
388  }
389  DTSegmentCand::AssPointCont pointsSet;
390 
391  // for (vector<DTSegmentCand::AssPoint>::const_iterator point=points.begin();
392  // point!=points.end(); ++point)
393  pointsSet.insert(points.begin(), points.end());
394 
395  if (debug) {
396  cout << "The Hits are" << endl;
397  copy(pointsSet.begin(), pointsSet.end(), ostream_iterator<DTSegmentCand::AssPoint>(std::cout));
398  cout << "----" << endl;
399  }
400 
401  DTSegmentCand* newCand = new DTSegmentCand(pointsSet, sl);
402  if (theUpdator->fit(newCand, false, false))
403  candidates.push_back(newCand);
404  else
405  delete newCand; // bad seg, too few hits
406  }
407 }
408 
410  for (vector<DTSegmentCand*>::iterator cand = cands.begin(); cand != cands.end(); ++cand)
411  if (*(*cand) == *seg)
412  return false;
413  return true;
414 }
std::pair< std::shared_ptr< DTHitPairForFit >, DTEnums::DTCellSide > AssPoint
Definition: DTSegmentCand.h:36
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
~DTCombinatorialPatternReco() override
Destructor.
std::vector< std::shared_ptr< DTHitPairForFit > > initHits(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
edm::OwnVector< DTSLRecSegment2D > reconstruct(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits) override
this function is called in the producer
DTCellSide
Which side of the DT cell.
Definition: DTEnums.h:15
std::vector< DTSegmentCand * > clean(const std::vector< DTSegmentCand *> &inputCands) const
do the cleaning
std::vector< std::shared_ptr< DTHitPairForFit > > hitCont
void update(DTRecSegment4D *seg, const bool calcT0, bool allow3par) const
recompute hits position and refit the segment4D
virtual bool good() const
void setES(const edm::EventSetup &setup)
set the setup
bool fit(DTSegmentCand *seg, bool allow3par, const bool fitdebug) const
edm::ESHandle< DTGeometry > theDTGeometry
void buildPointsCollection(std::vector< DTSegmentCand::AssPoint > &points, std::deque< std::shared_ptr< DTHitPairForFit >> &pointsNoLR, std::vector< DTSegmentCand *> &candidates, const DTSuperLayer *sl)
Basic3DVector unit() const
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
DTCombinatorialPatternReco(const edm::ParameterSet &pset, edm::ConsumesCollector)
Constructor.
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > theDTGeometryToken
std::set< AssPoint, AssPointLessZ > AssPointCont
Definition: DTSegmentCand.h:38
std::vector< DTSegmentCand * > buildSegments(const DTSuperLayer *sl, const std::vector< std::shared_ptr< DTHitPairForFit >> &hits)
hitCont::const_iterator hitIter
bool checkDoubleCandidates(std::vector< DTSegmentCand *> &segs, DTSegmentCand *seg)
DTSegmentCand * buildBestSegment(std::vector< DTSegmentCand::AssPoint > &assHits, const DTSuperLayer *sl)
std::vector< DTSegmentCand::AssPoint > findCompatibleHits(const LocalPoint &pos, const LocalVector &dir, const std::vector< std::shared_ptr< DTHitPairForFit >> &hits)
void setES(const edm::EventSetup &setup) override
DTSuperLayerId id() const
Return the DetId of this SL.
Definition: DTSuperLayer.cc:34
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72