CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTCombinatorialExtendedPatternReco.cc
Go to the documentation of this file.
1 
9 /* This Class Header */
11 
12 /* Collaborating Class Header */
27 
28 
29 /* C++ Headers */
30 #include <iterator>
31 using namespace std;
33 
34 /* ====================================================================== */
35 
38 DTRecSegment2DBaseAlgo(pset), theAlgoName("DTCombinatorialExtendedPatternReco")
39 {
40  theMaxAllowedHits = pset.getParameter<unsigned int>("MaxAllowedHits"); // 100
41  theAlphaMaxTheta = pset.getParameter<double>("AlphaMaxTheta");// 0.1 ;
42  theAlphaMaxPhi = pset.getParameter<double>("AlphaMaxPhi");// 1.0 ;
43  debug = pset.getUntrackedParameter<bool>("debug"); //true;
44  theUpdator = new DTSegmentUpdator(pset);
45  theCleaner = new DTSegmentCleaner(pset);
46  string theHitAlgoName = pset.getParameter<string>("recAlgo");
47  usePairs = !(theHitAlgoName=="DTNoDriftAlgo");
48 }
49 
52 }
53 
54 /* Operations */
57  const std::vector<DTRecHit1DPair>& pairs){
58 
59  if(debug) cout << "DTCombinatorialExtendedPatternReco::reconstruct" << endl;
60  theTriedPattern.clear();
62  vector<DTHitPairForFit*> hitsForFit = initHits(sl, pairs);
63 
64  vector<DTSegmentCand*> candidates = buildSegments(sl, hitsForFit);
65 
66  vector<DTSegmentCand*>::const_iterator cand=candidates.begin();
67  while (cand<candidates.end()) {
68 
69  DTSLRecSegment2D *segment = (**cand);
70 
71  theUpdator->update(segment);
72 
73  result.push_back(segment);
74 
75  if (debug) {
76  cout<<"Reconstructed 2D extended segments "<< result.back() <<endl;
77  }
78 
79  delete *(cand++); // delete the candidate!
80  }
81 
82  for (vector<DTHitPairForFit*>::iterator it = hitsForFit.begin(), ed = hitsForFit.end();
83  it != ed; ++it) delete *it;
84 
85  return result;
86 }
87 
89  // Get the DT Geometry
90  setup.get<MuonGeometryRecord>().get(theDTGeometry);
91  theUpdator->setES(setup);
92 }
93 
94 void DTCombinatorialExtendedPatternReco::setClusters(vector<DTSLRecCluster> clusters) {
95  theClusters = clusters;
96 }
97 
98 vector<DTHitPairForFit*>
100  const std::vector<DTRecHit1DPair>& hits){
101 
102  vector<DTHitPairForFit*> result;
103  for (vector<DTRecHit1DPair>::const_iterator hit=hits.begin();
104  hit!=hits.end(); ++hit) {
105  result.push_back(new DTHitPairForFit(*hit, *sl, theDTGeometry));
106  }
107  return result;
108 }
109 
110 vector<DTSegmentCand*>
112  const std::vector<DTHitPairForFit*>& hits){
113 
114  typedef vector<DTHitPairForFit*> hitCont;
115  typedef hitCont::const_iterator hitIter;
116  vector<DTSegmentCand*> result;
117 
118  if(debug) {
119  cout << "DTCombinatorialExtendedPatternReco::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  // hits must nor in the same nor in adiacent layers
143  if ( fabs((*lastHit)->id().layerId()-(*firstHit)->id().layerId())<=1 ) continue;
144  if(debug) {
145  cout << "Selected these two hits pair " << endl;
146  cout << "First " << *(*firstHit) << " Layer Id: " << (*firstHit)->id().layerId() << endl;
147  cout << "Last " << *(*lastHit) << " Layer Id: " << (*lastHit)->id().layerId() << endl;
148  }
149 
150  GlobalPoint IP;
151  float DAlphaMax;
152  if ((sl->id()).superlayer()==2) // Theta SL
153  DAlphaMax=theAlphaMaxTheta;
154  else // Phi SL
155  DAlphaMax=theAlphaMaxPhi;
156 
158  for (int firstLR=0; firstLR<2; ++firstLR) {
159  for (int lastLR=0; lastLR<2; ++lastLR) {
160  // TODO move the global transformation in the DTHitPairForFit class
161  // when it will be moved I will able to remove the sl from the input parameter
162  GlobalPoint gposFirst=sl->toGlobal( (*firstHit)->localPosition(codes[firstLR]) );
163  GlobalPoint gposLast= sl->toGlobal( (*lastHit)->localPosition(codes[lastLR]) );
164 
165  GlobalVector gvec=gposLast-gposFirst;
166  GlobalVector gvecIP=gposLast-IP;
167 
168  // difference in angle measured
169  float DAlpha=fabs(gvec.theta()-gvecIP.theta());
170 
171  // cout << "DAlpha " << DAlpha << endl;
172  if (DAlpha<DAlphaMax) {
173 
174  // create a segment hypotesis
175  // I don't need a true segment, just direction and position
176  LocalPoint posIni = (*firstHit)->localPosition(codes[firstLR]);
177  LocalVector dirIni =
178  ((*lastHit)->localPosition(codes[lastLR])-posIni).unit();
179 
180  // search for other compatible hits, with or without the L/R solved
181  vector<AssPoint> assHits = findCompatibleHits(posIni, dirIni, hits);
182  if(debug)
183  cout << "compatible hits " << assHits.size() << endl;
184 
185  // here return an extended candidate (which _has_ the original
186  // segment)
187  DTSegmentExtendedCand* seg = buildBestSegment(assHits, sl);
188 
189  if (seg) {
190  if(debug)
191  cout << "segment " << *seg<< endl;
192 
193  // check if the chi2 and #hits are ok
194  if (!seg->good()) { // good is reimplmented in extended segment
195  delete seg;
196  } else {
197 
198  // remove duplicated segments
199  if (checkDoubleCandidates(result,seg)) {
200  // add to the vector of hypotesis
201  // still work with extended segments
202  result.push_back(seg);
203  if(debug)
204  cout << "result is now " << result.size() << endl;
205  } else { // delete it!
206  delete seg;
207  if(debug)
208  cout << "already existing" << endl;
209  }
210  }
211  }
212  }
213  }
214  }
215  }
216  }
217  if (debug) {
218  for (vector<DTSegmentCand*>::const_iterator seg=result.begin();
219  seg!=result.end(); ++seg)
220  cout << *(*seg) << endl;
221  }
222 
223  // now I have a couple of segment hypotesis, should check for ghost
224  // still with extended candidates
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  // here, finally, I have to return the set of _original_ segments, not the
234  // extended ones.
235  return result;
236 }
237 
238 
239 vector<DTCombinatorialExtendedPatternReco::AssPoint>
241  const LocalVector& dirIni,
242  const vector<DTHitPairForFit*>& hits) {
243  if (debug) cout << "Pos: " << posIni << " Dir: "<< dirIni << endl;
244  vector<AssPoint> result;
245 
246  // counter to early-avoid double counting in hits pattern
247  vector<int> tried;
248  int nCompatibleHits=0;
249 
250  typedef vector<DTHitPairForFit*> hitCont;
251  typedef hitCont::const_iterator hitIter;
252  for (hitIter hit=hits.begin(); hit!=hits.end(); ++hit) {
253  pair<bool,bool> isCompatible = (*hit)->isCompatible(posIni, dirIni);
254  if (debug)
255  cout << "isCompatible " << isCompatible.first << " " <<
256  isCompatible.second << endl;
257 
258  // if only one of the two is compatible, then the LR is assigned,
259  // otherwise is undefined
260 
261  DTEnums::DTCellSide lrcode;
262  if (isCompatible.first && isCompatible.second) {
263  usePairs ? lrcode=DTEnums::undefLR : lrcode=DTEnums::Left ; // if not usePairs then only use single side
264  tried.push_back(3);
265  nCompatibleHits++;
266  }
267  else if (isCompatible.first) {
268  lrcode=DTEnums::Left;
269  tried.push_back(2);
270  nCompatibleHits++;
271  }
272  else if (isCompatible.second) {
273  lrcode=DTEnums::Right;
274  tried.push_back(1);
275  nCompatibleHits++;
276  }
277  else {
278  tried.push_back(0);
279  continue; // neither is compatible
280  }
281  result.push_back(AssPoint(*hit, lrcode));
282  }
283 
284 
285  // check if too few associated hits or pattern already tried
286  if ( nCompatibleHits < 3 || find(theTriedPattern.begin(), theTriedPattern.end(),tried) == theTriedPattern.end()) {
287  theTriedPattern.push_back(tried);
288  } else {
289  if (debug) {
290  vector<vector<int> >::const_iterator t=find(theTriedPattern.begin(),
291  theTriedPattern.end(),
292  tried);
293  cout << "Already tried";
294  copy((*t).begin(), (*t).end(), ostream_iterator<int>(std::cout));
295  cout << endl;
296  }
297  // empty the result vector
298  result.clear();
299  }
300  return result;
301 }
302 
305  const DTSuperLayer* sl) {
306  if (debug) cout << "DTCombinatorialExtendedPatternReco::buildBestSegment " <<
307  hits.size() << endl;
308  if (hits.size()<3) {
309  //cout << "buildBestSegment: hits " << hits.size()<< endl;
310  return 0; // a least 3 point
311  }
312 
313  // hits with defined LR
314  vector<AssPoint> points;
315 
316  // without: I store both L and R, a deque since I need front insertion and
317  // deletion
318  deque<DTHitPairForFit* > pointsNoLR;
319 
320  // first add only the hits with LR assigned
321  for (vector<AssPoint>::const_iterator hit=hits.begin();
322  hit!=hits.end(); ++hit) {
323  if ((*hit).second != DTEnums::undefLR) {
324  points.push_back(*hit);
325  } else { // then also for the undef'd one
326  pointsNoLR.push_back((*hit).first);
327  }
328  }
329 
330  if(debug) {
331  cout << "points " << points.size() << endl;
332  cout << "pointsNoLR " << pointsNoLR.size() << endl;
333  }
334 
335  // build all possible candidates using L/R ambiguity
336  vector<DTSegmentCand*> candidates ;
337 
338  buildPointsCollection(points, pointsNoLR, candidates, sl);
339 
340  // here I try to add the external clusters and build a set of "extended
341  // segment candidate
342  vector<DTSegmentExtendedCand*> extendedCands = extendCandidates(candidates,
343  sl);
344  if (debug) cout << "extended candidates " << extendedCands.size() << endl;
345 
346  // so now I have build a given number of segments, I should find the best one,
347  // by #hits and chi2.
348  vector<DTSegmentExtendedCand*>::const_iterator bestCandIter = extendedCands.end();
349  double minChi2=999999.;
350  unsigned int maxNumHits=0;
351  for (vector<DTSegmentExtendedCand*>::const_iterator iter=extendedCands.begin();
352  iter!=extendedCands.end(); ++iter) {
353  if ((*iter)->nHits()==maxNumHits && (*iter)->chi2()<minChi2) {
354  minChi2=(*iter)->chi2();
355  bestCandIter=iter;
356  } else if ((*iter)->nHits()>maxNumHits) {
357  maxNumHits=(*iter)->nHits();
358  minChi2=(*iter)->chi2();
359  bestCandIter=iter;
360  }
361  }
362 
363  // delete all candidates but the best one!
364  for (vector<DTSegmentExtendedCand*>::iterator iter=extendedCands.begin();
365  iter!=extendedCands.end(); ++iter)
366  if (iter!=bestCandIter) delete (*iter);
367 
368  // return the best candate if any
369  if (bestCandIter != extendedCands.end()) {
370  return (*bestCandIter);
371  }
372  return 0;
373 }
374 
375 void
377  deque<DTHitPairForFit*>& pointsNoLR,
378  vector<DTSegmentCand*>& candidates,
379  const DTSuperLayer* sl) {
380 
381  if(debug) {
382  cout << "DTCombinatorialExtendedPatternReco::buildPointsCollection " << endl;
383  cout << "points: " << points.size() << " NOLR: " << pointsNoLR.size()<< endl;
384  }
385  if (pointsNoLR.size()>0) { // still unassociated points!
386  DTHitPairForFit* unassHit = pointsNoLR.front();
387  // try with the right
388  if(debug)
389  cout << "Right hit" << endl;
390  points.push_back(AssPoint(unassHit, DTEnums::Right));
391  pointsNoLR.pop_front();
392  buildPointsCollection(points, pointsNoLR, candidates, sl);
393  pointsNoLR.push_front((unassHit));
394  points.pop_back();
395 
396  // try with the left
397  if(debug)
398  cout << "Left hit" << endl;
399  points.push_back(AssPoint(unassHit, DTEnums::Left));
400  pointsNoLR.pop_front();
401  buildPointsCollection(points, pointsNoLR, candidates, sl);
402  pointsNoLR.push_front((unassHit));
403  points.pop_back();
404  } else { // all associated
405 
406  if(debug) {
407  cout << "The Hits were" << endl;
408  copy(points.begin(), points.end(),
409  ostream_iterator<DTSegmentCand::AssPoint>(std::cout));
410  cout << "----" << endl;
411  cout << "All associated " << endl;
412  }
413  DTSegmentCand::AssPointCont pointsSet;
414 
415  // for (vector<AssPoint>::const_iterator point=points.begin();
416  // point!=points.end(); ++point)
417  pointsSet.insert(points.begin(),points.end());
418 
419  if(debug) {
420  cout << "The Hits are" << endl;
421  copy(pointsSet.begin(), pointsSet.end(),
422  ostream_iterator<DTSegmentCand::AssPoint>(std::cout));
423  cout << "----" << endl;
424  }
425 
426  DTSegmentCand* newCand = new DTSegmentCand(pointsSet,sl);
427  if (theUpdator->fit(newCand)) candidates.push_back(newCand);
428  else delete newCand; // bad seg, too few hits
429  }
430 }
431 
432 bool
434  DTSegmentCand* seg) {
435  for (vector<DTSegmentCand*>::iterator cand=cands.begin();
436  cand!=cands.end(); ++cand)
437  if (*(*cand)==*seg) return false;
438  return true;
439 }
440 
441 vector<DTSegmentExtendedCand*>
442 DTCombinatorialExtendedPatternReco::extendCandidates(vector<DTSegmentCand*>& candidates,
443  const DTSuperLayer* sl) {
444  if (debug) cout << "extendCandidates " << candidates.size() << endl;
445  vector<DTSegmentExtendedCand*> result;
446 
447  // in case of phi SL just return
448  if (sl->id().superLayer() != 2 ) {
449  for (vector<DTSegmentCand*>:: const_iterator cand=candidates.begin();
450  cand!=candidates.end(); ++cand) {
451  DTSegmentExtendedCand* extendedCand = new DTSegmentExtendedCand(*cand);
452  // and delete the original candidate
453  delete *cand;
454  result.push_back(extendedCand);
455  }
456  return result;
457  }
458 
459  // first I have to select the cluster which are compatible with the actual
460  // candidate, namely +/-1 sector/station/wheel
461  vector<DTSegmentExtendedCand::DTSLRecClusterForFit> clustersWithPos;
462  if (debug) cout << "AllClustersWithPos " << theClusters.size() << endl;
463  if(debug) cout << "SL: " << sl->id() << endl;
464  for (vector<DTSLRecCluster>::const_iterator clus=theClusters.begin();
465  clus!=theClusters.end(); ++clus) {
466  if(debug) cout << "CLUS: " << (*clus).superLayerId() << endl;
467  if ((*clus).superLayerId().superLayer()==2 && closeSL(sl->id(),(*clus).superLayerId())) {
468  // and then get their pos in the actual SL frame
469  const DTSuperLayer* clusSl =
470  theDTGeometry->superLayer((*clus).superLayerId());
471  LocalPoint pos=sl->toLocal(clusSl->toGlobal((*clus).localPosition()));
472  //LocalError err=sl->toLocal(clusSl->toGlobal((*clus).localPositionError()));
473  LocalError err=(*clus).localPositionError();
474  clustersWithPos.push_back(DTSegmentExtendedCand::DTSLRecClusterForFit(*clus, pos, err) );
475  }
476  }
477  if (debug) cout << "closeClustersWithPos " << clustersWithPos.size() << endl;
478 
479  for (vector<DTSegmentCand*>:: const_iterator cand=candidates.begin();
480  cand!=candidates.end(); ++cand) {
481  // create an extended candidate
482  DTSegmentExtendedCand* extendedCand = new DTSegmentExtendedCand(*cand);
483  // and delete the original candidate
484  delete *cand;
485  // do this only for theta SL
486  if (extendedCand->superLayer()->id().superLayer() == 2 ) {
487  // first check compatibility between cand and clusForFit
488  for (vector<DTSegmentExtendedCand::DTSLRecClusterForFit>::const_iterator
489  exClus=clustersWithPos.begin(); exClus!=clustersWithPos.end(); ++exClus) {
490  if (extendedCand->isCompatible(*exClus)) {
491  if (debug) cout << "is compatible " << endl;
492  // add compatible cluster
493  extendedCand->addClus(*exClus);
494  }
495  }
496  // fit the segment
497  if (debug) cout << "extended cands nHits: " << extendedCand->nHits() <<endl;
498  if (theUpdator->fit(extendedCand)) {
499  // add to result
500  result.push_back(extendedCand);
501  } else {
502  cout << "Bad fit" << endl;
503  delete extendedCand;
504  }
505  } else { // Phi SuperLayer
506  result.push_back(extendedCand);
507  }
508  }
509 
510  return result;
511 }
512 
514  const DTSuperLayerId& id2) {
515  if (id1==id2) return false;
516  if (abs(id1.wheel()-id2.wheel())>1 ) return false;
517  // take into account also sector 13 and 14
518  int sec1 = ( id1.sector()==13 ) ? 4: id1.sector();
519  sec1=(sec1==14)? 10: sec1;
520  int sec2 = ( id2.sector()==13 ) ? 4: id2.sector();
521  sec2=(sec2==14)? 10: sec2;
522  // take into account also sector 1/12
523  if (abs(sec1-sec2)>1 && abs(sec1-sec2)!=11 ) return false;
524  //if (abs(id1.station()-id2.station())>1 ) return false;
525  return true;
526 }
T getParameter(std::string const &) const
std::vector< DTHitPairForFit * > initHits(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
T getUntrackedParameter(std::string const &, T const &) const
std::vector< AssPoint > findCompatibleHits(const LocalPoint &pos, const LocalVector &dir, const std::vector< DTHitPairForFit * > &hits)
reference back()
Definition: OwnVector.h:320
std::vector< DTSegmentCand * > buildSegments(const DTSuperLayer *sl, const std::vector< DTHitPairForFit * > &hits)
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
virtual edm::OwnVector< DTSLRecSegment2D > reconstruct(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
this function is called in the producer
std::pair< DTHitPairForFit *, DTEnums::DTCellSide > AssPoint
#define abs(x)
Definition: mlp_lapack.h:159
DTCellSide
Which side of the DT cell.
Definition: DTEnums.h:17
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:62
DTCombinatorialExtendedPatternReco(const edm::ParameterSet &pset)
Constructor.
bool checkDoubleCandidates(std::vector< DTSegmentCand * > &segs, DTSegmentCand *seg)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
virtual void setES(const edm::EventSetup &setup)
virtual bool good() const
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
DTSegmentExtendedCand * buildBestSegment(std::vector< AssPoint > &assHits, const DTSuperLayer *sl)
DTSuperLayerId id() const
Return the DetId of this SL.
Definition: DTSuperLayer.cc:38
void push_back(D *&d)
Definition: OwnVector.h:273
void buildPointsCollection(std::vector< AssPoint > &points, std::deque< DTHitPairForFit * > &pointsNoLR, std::vector< DTSegmentCand * > &candidates, const DTSuperLayer *sl)
void setES(const edm::EventSetup &setup)
set the setup
string unit
Definition: csvLumiCalc.py:46
std::vector< std::vector< int > > theTriedPattern
tuple result
Definition: query.py:137
tuple IP
Definition: listHistos.py:76
int superLayer() const
Return the superlayer number.
bool fit(DTSegmentCand *seg) const
void addClus(const DTSegmentExtendedCand::DTSLRecClusterForFit &clus)
std::vector< DTSegmentExtendedCand * > extendCandidates(std::vector< DTSegmentCand * > &candidates, const DTSuperLayer *sl)
void setClusters(std::vector< DTSLRecCluster > clusters)
const DTSuperLayer * superLayer() const
the super layer on which relies
Definition: DTSegmentCand.h:77
bool isCompatible(const DTSegmentExtendedCand::DTSLRecClusterForFit &clus)
virtual unsigned int nHits() const
std::set< AssPoint, AssPointLessZ > AssPointCont
Definition: DTSegmentCand.h:41
const T & get() const
Definition: EventSetup.h:55
bool closeSL(const DTSuperLayerId &id1, const DTSuperLayerId &id2)
std::vector< DTSegmentCand * > clean(std::vector< DTSegmentCand * > inputCands) const
do the cleaning
int sector() const
Definition: DTChamberId.h:63
tuple cout
Definition: gather_cfg.py:121
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")