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