CMS 3D CMS Logo

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