CMS 3D CMS Logo

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