CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTMeantimerPatternReco.cc
Go to the documentation of this file.
1 
10 /* This Class Header */
12 
13 /* Collaborating Class Header */
27 
28 /* C++ Headers */
29 #include <iterator>
30 using namespace std;
32 
33 /* ====================================================================== */
34 
35 
38 DTRecSegment2DBaseAlgo(pset), theAlgoName("DTMeantimerPatternReco")
39 {
40 
41  theMaxAllowedHits = pset.getParameter<unsigned int>("MaxAllowedHits"); // 100
42  theAlphaMaxTheta = pset.getParameter<double>("AlphaMaxTheta");// 0.1 ;
43  theAlphaMaxPhi = pset.getParameter<double>("AlphaMaxPhi");// 1.0 ;
44  theMaxT0 = pset.getParameter<double>("MaxT0");
45  theMinT0 = pset.getParameter<double>("MinT0");
46  theMaxChi2 = pset.getParameter<double>("MaxChi2");// 8.0 ;
47  debug = pset.getUntrackedParameter<bool>("debug");
48  theUpdator = new DTSegmentUpdator(pset);
49  theCleaner = new DTSegmentCleaner(pset);
50 }
51 
54  delete theUpdator;
55  delete theCleaner;
56 }
57 
58 /* Operations */
61  const std::vector<DTRecHit1DPair>& pairs){
62 
64  vector<DTHitPairForFit*> hitsForFit = initHits(sl, pairs);
65 
66  vector<DTSegmentCand*> candidates = buildSegments(sl, hitsForFit);
67 
68  vector<DTSegmentCand*>::const_iterator cand=candidates.begin();
69  while (cand<candidates.end()) {
70 
71  DTSLRecSegment2D *segment = (**cand);
72 
73  theUpdator->update(segment);
74 
75  if (debug)
76  cout<<"Reconstructed 2D segments "<<*segment<<endl;
77  result.push_back(segment);
78 
79  delete *(cand++); // delete the candidate!
80  }
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 
107  typedef vector<DTHitPairForFit*> hitCont;
108  typedef hitCont::const_iterator hitIter;
109  vector<DTSegmentCand*> result;
111 
112  if(debug) {
113  cout << "buildSegments: " << sl->id() << " nHits " << hits.size() << endl;
114  for (hitIter hit=hits.begin(); hit!=hits.end(); ++hit) cout << **hit<< " wire: " << (*hit)->id() << " DigiTime: " << (*hit)->digiTime() << endl;
115  }
116 
117  // 10-Mar-2004 SL
118  // put a protection against heavily populated chambers, for which the segment
119  // building could lead to infinite memory usage...
120  if (hits.size() > theMaxAllowedHits ) {
121  if(debug) {
122  cout << "Warning: this SuperLayer " << sl->id() << " has too many hits : "
123  << hits.size() << " max allowed is " << theMaxAllowedHits << endl;
124  cout << "Skipping segment reconstruction... " << endl;
125  }
126  return result;
127  }
128 
129  GlobalPoint IP;
130  float DAlphaMax;
131  if ((sl->id()).superlayer()==2) // Theta SL
132  DAlphaMax=theAlphaMaxTheta;
133  else // Phi SL
134  DAlphaMax=theAlphaMaxPhi;
135 
136  vector<AssPoint> usedHits;
137 
138  // get two hits in different layers and see if there are other hits
139  // compatible with them
140  for (hitCont::const_iterator firstHit=hits.begin(); firstHit!=hits.end();
141  ++firstHit) {
142  for (hitCont::const_reverse_iterator lastHit=hits.rbegin();
143  (*lastHit)!=(*firstHit); ++lastHit) {
144 
145  // a geometrical sensibility cut for the two hits
146  if (!geometryFilter((*firstHit)->id(),(*lastHit)->id())) continue;
147 
148  // create a set of hits for the fit (only the hits between the two selected ones)
149  hitCont hitsForFit;
150  for (hitCont::const_iterator tmpHit=firstHit+1; (*tmpHit)!=(*lastHit); tmpHit++)
151  if ((geometryFilter((*tmpHit)->id(),(*lastHit)->id()))
152  && (geometryFilter((*tmpHit)->id(),(*firstHit)->id()))) hitsForFit.push_back(*tmpHit);
153 
154  for (int firstLR=0; firstLR<2; ++firstLR) {
155  for (int lastLR=0; lastLR<2; ++lastLR) {
156 
157  // TODO move the global transformation in the DTHitPairForFit class
158  // when it will be moved I will able to remove the sl from the input parameter
159  GlobalPoint gposFirst=sl->toGlobal( (*firstHit)->localPosition(codes[firstLR]) );
160  GlobalPoint gposLast= sl->toGlobal( (*lastHit)->localPosition(codes[lastLR]) );
161  GlobalVector gvec=gposLast-gposFirst;
162  GlobalVector gvecIP=gposLast-IP;
163 
164  // difference in angle measured
165  float DAlpha=fabs(gvec.theta()-gvecIP.theta());
166  if (DAlpha>DAlphaMax) continue;
167 
168  if(debug) {
169  cout << "Selected hit pair:" << endl;
170  cout << " First " << *(*firstHit) << " Layer Id: " << (*firstHit)->id().layerId() << " Side: " << firstLR << " DigiTime: " << (*firstHit)->digiTime() << endl;
171  cout << " Last " << *(*lastHit) << " Layer Id: " << (*lastHit)->id().layerId() << " Side: " << lastLR << " DigiTime: " << (*lastHit)->digiTime() << endl;
172  }
173 
174  vector<AssPoint> assHits;
175  // create a candidate hit list
176  assHits.push_back(AssPoint(*firstHit,codes[firstLR]));
177  assHits.push_back(AssPoint(*lastHit,codes[lastLR]));
178 
179  // run hit adding/segment building
180  maxfound = 3;
181  addHits(sl,assHits,hitsForFit,result, usedHits);
182  }
183  }
184  }
185  }
186 
187  // now I have a couple of segment hypotheses, should check for ghost
188  if (debug) {
189  cout << "Result (before cleaning): " << result.size() << endl;
190  for (vector<DTSegmentCand*>::const_iterator seg=result.begin(); seg!=result.end(); ++seg)
191  cout << *(*seg) << endl;
192  }
193  result = theCleaner->clean(result);
194  if (debug) {
195  cout << "Result (after cleaning): " << result.size() << endl;
196  for (vector<DTSegmentCand*>::const_iterator seg=result.begin();
197  seg!=result.end(); ++seg)
198  cout << *(*seg) << endl;
199  }
200 
201  return result;
202 }
203 
204 
205 
206 void
207 DTMeantimerPatternReco::addHits(const DTSuperLayer* sl, vector<AssPoint>& assHits, const vector<DTHitPairForFit*>& hits,
208  vector<DTSegmentCand*> &result, vector<AssPoint>& usedHits) {
209 
210  typedef vector<DTHitPairForFit*> hitCont;
211  double chi2l,chi2r,t0_corrl,t0_corrr;
212  bool foundSomething = false;
213 
214  if (debug)
215  cout << "DTMeantimerPatternReco::addHit " << endl << " Picked " << assHits.size() << " hits, " << hits.size() << " left." << endl;
216 
217  if (assHits.size()+hits.size()<maxfound) return;
218 
219  // loop over the remaining hits
220  for (hitCont::const_iterator hit=hits.begin(); hit!=hits.end(); ++hit) {
221  if (debug)
222  cout << " Trying B: " << **hit<< " wire: " << (*hit)->id() << endl;
223 
224  assHits.push_back(AssPoint(*hit, DTEnums::Left));
225  bool left_ok=fitWithT0(assHits, chi2l, t0_corrl,0);
226  assHits.pop_back();
227 
228  assHits.push_back(AssPoint(*hit, DTEnums::Right));
229  bool right_ok=fitWithT0(assHits, chi2r, t0_corrr,0);
230  assHits.pop_back();
231 
232  if (debug) {
233  int nHits=assHits.size()+1;
234  cout << " Left: t0_corr = " << t0_corrl << "ns chi2/nHits = " << chi2l << "/" << nHits << " ok: " << left_ok << endl;
235  cout << " Right: t0_corr = " << t0_corrr << "ns chi2/nHits = " << chi2r << "/" << nHits << " ok: " << right_ok << endl;
236  }
237 
238  if (!left_ok && !right_ok) continue;
239 
240  foundSomething = true;
241 
242  // prepare the hit set for the next search, start from the other side
243  hitCont hitsForFit;
244  for (hitCont::const_iterator tmpHit=hit+1; tmpHit!=hits.end(); tmpHit++)
245  if (geometryFilter((*tmpHit)->id(),(*hit)->id())) hitsForFit.push_back(*tmpHit);
246 
247  reverse(hitsForFit.begin(),hitsForFit.end());
248 
249  // choose only one - left or right
250  if (assHits.size()>3 && left_ok && right_ok) {
251  if (chi2l<chi2r-0.1) right_ok=false; else
252  if (chi2r<chi2l-0.1) left_ok=false;
253  }
254  if (left_ok) {
255  assHits.push_back(AssPoint(*hit, DTEnums::Left));
256  addHits(sl,assHits,hitsForFit,result,usedHits);
257  assHits.pop_back();
258  }
259  if (right_ok) {
260  assHits.push_back(AssPoint(*hit, DTEnums::Right));
261  addHits(sl,assHits,hitsForFit,result,usedHits);
262  assHits.pop_back();
263  }
264  }
265 
266  if (foundSomething) return;
267 
268  // If we already have a segment with more hits from this hit pair - don't save this one.
269  if (assHits.size()<maxfound) return;
270 
271  // Check if semgent Ok, calculate chi2
272  if (!fitWithT0(assHits, chi2l, t0_corrl,debug)) return;
273 
274  // If no more iterations - store the current segment
275 
276  DTSegmentCand::AssPointCont pointsSet;
277  pointsSet.insert(assHits.begin(),assHits.end());
278  DTSegmentCand* seg = new DTSegmentCand(pointsSet,sl);
279  theUpdator->fit(seg);
280 
281  if (seg) {
282  for (vector<AssPoint>::const_iterator hiti = assHits.begin()+1; hiti != assHits.end()-1; hiti++)
283  usedHits.push_back(*hiti);
284 
285  if (assHits.size()>maxfound) maxfound = assHits.size();
286  seg->setChi2(chi2l); // we need to set the chi^2 so that the cleaner can pick the best segments
287  if (debug) cout << "Segment built: " << *seg<< endl;
288  if (checkDoubleCandidates(result,seg)) {
289  result.push_back(seg);
290  if (debug) cout << " Result is now " << result.size() << endl;
291  } else {
292  if (debug) cout << " Exists - skipping" << endl;
293  delete seg;
294  }
295  }
296 }
297 
298 
299 bool
301 {
302 // return true;
303 
304  const int layerLowerCut[4]={0,-1,-2,-2};
305  const int layerUpperCut[4]={0, 2, 2, 3};
306 // const int layerLowerCut[4]={0,-2,-4,-5};
307 // const int layerUpperCut[4]={0, 3, 4, 6};
308 
309  // deal only with hits that are in the same SL
310  if (first.layerId().superlayerId().superLayer()!=second.layerId().superlayerId().superLayer())
311  return true;
312 
313  int deltaLayer=abs(first.layerId().layer()-second.layerId().layer());
314 
315  // drop hits in the same layer
316  if (!deltaLayer) return false;
317 
318  // protection against unexpected layer numbering
319  if (deltaLayer>3) {
320  cout << "*** WARNING! DT Layer numbers differ by more than 3! for hits: " << endl;
321  cout << " " << first << endl;
322  cout << " " << second << endl;
323  return false;
324  }
325 
326  // accept only hits in cells "not too far away"
327  int deltaWire=first.wire()-second.wire();
328  if (second.layerId().layer()%2==0) deltaWire=-deltaWire; // yet another trick to get it right...
329  if ((deltaWire<=layerLowerCut[deltaLayer]) || (deltaWire>=layerUpperCut[deltaLayer])) return false;
330 
331  return true;
332 }
333 
334 
335 bool
336 DTMeantimerPatternReco::fitWithT0(const vector<AssPoint> &assHits, double &chi2, double &t0_corr, const bool fitdebug)
337 {
338  typedef vector < pair<double,double> > hitCoord;
339  double a,b,coordError,x,y;
340  double sx=0,sy=0,sxy=0,sxx=0,ssx=0,ssy=0,s=0,ss=0;
341  int leftHits=0,rightHits=0;
342 
343  if (assHits.size()<3) return false;
344 
345  // I'm assuming the single hit error is the same for all hits...
346  coordError=((*(assHits.begin())).first)->localPositionError().xx();
347 
348  for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
349  if (coordError!=((*(assHits.begin())).first)->localPositionError().xx())
350  cout << " DTMeantimerPatternReco: Warning! Hit errors different within one segment!" << endl;
351 
352  x=((*assHit).first)->localPosition((*assHit).second).z();
353  y=((*assHit).first)->localPosition((*assHit).second).x();
354 
355  sx+=x;
356  sy+=y;
357  sxy+=x*y;
358  sxx+=x*x;
359  s++;
360 
361  if ((*assHit).second==DTEnums::Left) {
362  leftHits++;
363  ssx+=x;
364  ssy+=y;
365  ss++;
366  } else {
367  rightHits++;
368  ssx-=x;
369  ssy-=y;
370  ss--;
371  }
372  }
373 
374  if (fitdebug && debug)
375  cout << " DTMeantimerPatternReco::fitWithT0 Left hits: " << leftHits << " Right hits: " << rightHits << endl;
376 
377  if (leftHits && rightHits) {
378 
379  double delta = ss*ss*sxx+s*sx*sx+s*ssx*ssx-s*s*sxx-2*ss*sx*ssx;
380  t0_corr=0.;
381 
382  if (delta) {
383  a=(ssy*s*ssx+sxy*ss*ss+sy*sx*s-sy*ss*ssx-ssy*sx*ss-sxy*s*s)/delta;
384  b=(ssx*sy*ssx+sxx*ssy*ss+sx*sxy*s-sxx*sy*s-ssx*sxy*ss-sx*ssy*ssx)/delta;
385  t0_corr=(ssx*s*sxy+sxx*ss*sy+sx*sx*ssy-sxx*s*ssy-sx*ss*sxy-ssx*sx*sy)/delta;
386  } else return false;
387  } else {
388  double d = s*sxx - sx*sx;
389  b = (sxx*sy- sx*sxy)/ d;
390  a = (s*sxy - sx*sy) / d;
391  t0_corr=0;
392  }
393 
394 // Calculate the chi^2 of the hits AFTER t0 correction
395 
396  double chi,chi2_not0;
397  chi2=0;
398  chi2_not0=0;
399 
400  for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
401  x=((*assHit).first)->localPosition((*assHit).second).z();
402  y=((*assHit).first)->localPosition((*assHit).second).x();
403 
404  chi=y-a*x-b;
405  chi2_not0+=chi*chi/coordError;
406 
407  if ((*assHit).second==DTEnums::Left) chi-=t0_corr;
408  else chi+=t0_corr;
409 
410  chi2+=chi*chi/coordError;
411  }
412 
413  // For 3-hit segments ignore timing information
414  if (assHits.size()<4) {
415  chi2=chi2_not0;
416 // if (chi2<theMaxChi2) return true;
417  if (chi2<200.) return true;
418  else return false;
419  }
420 
421  t0_corr/=-0.00543; // convert drift distance to time
422 
423  if (debug && fitdebug)
424  {
425  cout << " t0_corr = " << t0_corr << "ns chi2/nHits = " << chi2 << "/" << assHits.size() << endl;
426 
427  if (((chi2/(assHits.size()-2)<theMaxChi2) && (t0_corr<theMaxT0) && (t0_corr>theMinT0)) || (assHits.size()==4)) {
428  cout << " a = " << a << " b = " << b << endl;
429  for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
430  x=((*assHit).first)->localPosition((*assHit).second).z();
431  y=((*assHit).first)->localPosition((*assHit).second).x();
432  cout << " z= " << x << " x= " << y;
433  if ((*assHit).second==DTEnums::Left) cout << " x_corr= " << y+t0_corr*0.00543;
434  else cout << " x_corr= " << y-t0_corr*0.00543;
435  cout << " seg= " << a*x+b << endl;
436  }
437  }
438  }
439 
440  if ((chi2/(assHits.size()-2)<theMaxChi2) && (t0_corr<theMaxT0) && (t0_corr>theMinT0)) return true;
441  else return false;
442 }
443 
444 void
445 DTMeantimerPatternReco::rawFit(double &a, double &b, const vector< pair<double,double> > &hits) {
446 
447  double s=0,sx=0,sy=0,x,y;
448  double sxx=0,sxy=0;
449 
450  a=b=0;
451  if (hits.size()==0) return;
452 
453  if (hits.size()==1) {
454  b=(*(hits.begin())).second;
455  } else {
456  for (unsigned int i = 0; i != hits.size(); i++) {
457  x=hits[i].first;
458  y=hits[i].second;
459  sy += y;
460  sxy+= x*y;
461  s += 1.;
462  sx += x;
463  sxx += x*x;
464  }
465  // protect against a vertical line???
466 
467  double d = s*sxx - sx*sx;
468  b = (sxx*sy- sx*sxy)/ d;
469  a = (s*sxy - sx*sy) / d;
470  }
471 }
472 
473 bool
474 DTMeantimerPatternReco::checkDoubleCandidates(vector<DTSegmentCand*>& cands,
475  DTSegmentCand* seg) {
476  for (vector<DTSegmentCand*>::iterator cand=cands.begin();
477  cand!=cands.end(); ++cand) {
478  if (*(*cand)==*seg) return false;
479  if (((*cand)->nHits()>=seg->nHits()) && ((*cand)->chi2ndof()<seg->chi2ndof()))
480  if ((*cand)->nSharedHitPairs(*seg)>int(seg->nHits()-2)) return false;
481  }
482  return true;
483 }
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
std::pair< DTHitPairForFit *, DTEnums::DTCellSide > AssPoint
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void update(DTRecSegment4D *seg, const bool calcT0=false) const
recompute hits position and refit the segment4D
virtual double chi2ndof() const
the chi2/NDOF of the fit
Definition: DTSegmentCand.h:68
bool checkDoubleCandidates(std::vector< DTSegmentCand * > &segs, DTSegmentCand *seg)
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
#define abs(x)
Definition: mlp_lapack.h:159
DTCellSide
Which side of the DT cell.
Definition: DTEnums.h:17
virtual void setChi2(double &chi2)
set chi2
virtual ~DTMeantimerPatternReco()
Destructor.
int layer() const
Return the layer number.
Definition: DTLayerId.h:55
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:61
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
U second(std::pair< T, U > const &p)
DTSuperLayerId id() const
Return the DetId of this SL.
Definition: DTSuperLayer.cc:38
void push_back(D *&d)
Definition: OwnVector.h:273
virtual void setES(const edm::EventSetup &setup)
void setES(const edm::EventSetup &setup)
set the setup
std::vector< DTHitPairForFit * > initHits(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
tuple result
Definition: query.py:137
tuple IP
Definition: listHistos.py:76
int superLayer() const
Return the superlayer number.
bool fit(DTSegmentCand *seg) const
bool fitWithT0(const std::vector< AssPoint > &assHits, double &chi2, double &t0_corr, const bool fitdebug)
bool first
Definition: L1TdeRCT.cc:94
unsigned int id
int wire() const
Return the wire number.
Definition: DTWireId.h:58
std::set< AssPoint, AssPointLessZ > AssPointCont
Definition: DTSegmentCand.h:41
void rawFit(double &a, double &b, const std::vector< std::pair< double, double > > &hits)
virtual unsigned int nHits() const
Definition: DTSegmentCand.h:62
std::vector< DTSegmentCand * > buildSegments(const DTSuperLayer *sl, const std::vector< DTHitPairForFit * > &hits)
const T & get() const
Definition: EventSetup.h:55
double b
Definition: hdecay.h:120
bool geometryFilter(const DTWireId first, const DTWireId second) const
virtual edm::OwnVector< DTSLRecSegment2D > reconstruct(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
this function is called in the producer
std::vector< DTSegmentCand * > clean(std::vector< DTSegmentCand * > inputCands) const
do the cleaning
void addHits(const DTSuperLayer *sl, std::vector< AssPoint > &assHits, const std::vector< DTHitPairForFit * > &hits, std::vector< DTSegmentCand * > &result, std::vector< AssPoint > &usedHits)
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:64
DTMeantimerPatternReco(const edm::ParameterSet &pset)
Constructor.
double a
Definition: hdecay.h:121
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
edm::ESHandle< DTGeometry > theDTGeometry