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 
8 /* This Class Header */
10 
11 /* Collaborating Class Header */
26 
27 /* C++ Headers */
28 #include <iterator>
29 using namespace std;
31 
32 /* ====================================================================== */
33 
34  typedef std::vector<std::shared_ptr<DTHitPairForFit>> hitCont;
35  typedef hitCont::const_iterator hitIter;
36 
37 
41  theFitter(new DTLinearFit()),
42  theAlgoName("DTMeantimerPatternReco")
43 {
44  theMaxAllowedHits = pset.getParameter<unsigned int>("MaxAllowedHits"); // 100
45  theAlphaMaxTheta = pset.getParameter<double>("AlphaMaxTheta");// 0.1 ;
46  theAlphaMaxPhi = pset.getParameter<double>("AlphaMaxPhi");// 1.0 ;
47  theMaxChi2 = pset.getParameter<double>("MaxChi2");// 8.0 ;
48  debug = pset.getUntrackedParameter<bool>("debug");
49  theUpdator = new DTSegmentUpdator(pset);
50  theCleaner = new DTSegmentCleaner(pset);
51 }
52 
55  delete theUpdator;
56  delete theCleaner;
57  delete theFitter;
58 }
59 
60 /* Operations */
63  const std::vector<DTRecHit1DPair>& pairs){
64 
66  std::vector<std::shared_ptr<DTHitPairForFit>> hitsForFit = initHits(sl, pairs);
67 
68  vector<DTSegmentCand*> candidates = buildSegments(sl, hitsForFit);
69 
70  vector<DTSegmentCand*>::const_iterator cand=candidates.begin();
71  while (cand<candidates.end()) {
72 
73  DTSLRecSegment2D *segment = (**cand);
74  theUpdator->update(segment,1);
75 
76  if (debug)
77  cout<<"Reconstructed 2D segments "<<*segment<<endl;
78  result.push_back(segment);
79 
80  delete *(cand++); // delete the candidate!
81  }
82 
83  return result;
84 }
85 
87  // Get the DT Geometry
88  setup.get<MuonGeometryRecord>().get(theDTGeometry);
89  theUpdator->setES(setup);
90 }
91 
92 vector<std::shared_ptr<DTHitPairForFit>>
94  const std::vector<DTRecHit1DPair>& hits){
95 
97  for (vector<DTRecHit1DPair>::const_iterator hit=hits.begin();
98  hit!=hits.end(); ++hit) {
99  result.push_back(std::make_shared<DTHitPairForFit>(*hit, *sl, theDTGeometry));
100  }
101  return result;
102 }
103 
104 vector<DTSegmentCand*>
106  const std::vector<std::shared_ptr<DTHitPairForFit>>& hits){
107 
108  vector<DTSegmentCand*> result;
110 
111  if(debug) {
112  cout << "buildSegments: " << sl->id() << " nHits " << hits.size() << endl;
113  for (hitIter hit=hits.begin(); hit!=hits.end(); ++hit) cout << **hit<< " wire: " << (*hit)->id() << " DigiTime: " << (*hit)->digiTime() << endl;
114  }
115 
116  if (hits.size() > theMaxAllowedHits ) {
117  if(debug) {
118  cout << "Warning: this SuperLayer " << sl->id() << " has too many hits : "
119  << hits.size() << " max allowed is " << theMaxAllowedHits << endl;
120  cout << "Skipping segment reconstruction... " << endl;
121  }
122  return result;
123  }
124 
125  GlobalPoint IP;
126  float DAlphaMax;
127  if ((sl->id()).superlayer()==2) // Theta SL
128  DAlphaMax=theAlphaMaxTheta;
129  else // Phi SL
130  DAlphaMax=theAlphaMaxPhi;
131 
132  // get two hits in different layers and see if there are other hits
133  // compatible with them
134  for (hitCont::const_iterator firstHit=hits.begin(); firstHit!=hits.end();
135  ++firstHit) {
136  for (hitCont::const_reverse_iterator lastHit=hits.rbegin();
137  (*lastHit)!=(*firstHit); ++lastHit) {
138 
139  // a geometrical sensibility cut for the two hits
140  if (!geometryFilter((*firstHit)->id(),(*lastHit)->id())) continue;
141 
142  // create a set of hits for the fit (only the hits between the two selected ones)
143  hitCont hitsForFit;
144  for (hitCont::const_iterator tmpHit=firstHit+1; (*tmpHit)!=(*lastHit); tmpHit++)
145  if ((geometryFilter((*tmpHit)->id(),(*lastHit)->id()))
146  && (geometryFilter((*tmpHit)->id(),(*firstHit)->id()))) hitsForFit.push_back(*tmpHit);
147 
148  for (int firstLR=0; firstLR<2; ++firstLR) {
149  for (int lastLR=0; lastLR<2; ++lastLR) {
150 
151  // TODO move the global transformation in the DTHitPairForFit class
152  // when it will be moved I will able to remove the sl from the input parameter
153  GlobalPoint gposFirst=sl->toGlobal( (*firstHit)->localPosition(codes[firstLR]) );
154  GlobalPoint gposLast= sl->toGlobal( (*lastHit)->localPosition(codes[lastLR]) );
155  GlobalVector gvec=gposLast-gposFirst;
156  GlobalVector gvecIP=gposLast-IP;
157 
158  // difference in angle measured
159  float DAlpha=fabs(gvec.theta()-gvecIP.theta());
160  if (DAlpha>DAlphaMax) continue;
161 
162 // if(debug) {
163 // cout << "Selected hit pair:" << endl;
164 // cout << " First " << *(*firstHit) << " Layer Id: " << (*firstHit)->id().layerId() << " Side: " << firstLR << " DigiTime: " << (*firstHit)->digiTime() << endl;
165 // cout << " Last " << *(*lastHit) << " Layer Id: " << (*lastHit)->id().layerId() << " Side: " << lastLR << " DigiTime: " << (*lastHit)->digiTime() << endl;
166 // }
167 
168  vector<DTSegmentCand::AssPoint> assHits;
169  // create a candidate hit list
170  assHits.push_back(DTSegmentCand::AssPoint(*firstHit,codes[firstLR]));
171  assHits.push_back(DTSegmentCand::AssPoint(*lastHit,codes[lastLR]));
172 
173  // run hit adding/segment building
174  maxfound = 3;
175  addHits(sl,assHits,hitsForFit,result);
176  }
177  }
178  }
179  }
180 
181  // now I have a couple of segment hypotheses, should check for ghosts
182  if (debug) {
183  cout << "Result (before cleaning): " << result.size() << endl;
184  for (vector<DTSegmentCand*>::const_iterator seg=result.begin(); seg!=result.end(); ++seg) cout << *(*seg) << endl;
185  }
186 
187  result = theCleaner->clean(result);
188 
189  if (debug) {
190  cout << "Result (after cleaning): " << result.size() << endl;
191  for (vector<DTSegmentCand*>::const_iterator seg=result.begin(); seg!=result.end(); ++seg) cout << *(*seg) << endl;
192  }
193 
194  return result;
195 }
196 
197 
198 
199 void
200 DTMeantimerPatternReco::addHits(const DTSuperLayer* sl, vector<DTSegmentCand::AssPoint>& assHits, const vector<std::shared_ptr<DTHitPairForFit>>& hits,
201  vector<DTSegmentCand*> &result) {
202 
203  double chi2l,chi2r,t0_corrl,t0_corrr;
204  bool foundSomething = false;
205 
206 // if (debug)
207 // cout << "DTMeantimerPatternReco::addHit " << endl << " Picked " << assHits.size() << " hits, " << hits.size() << " left." << endl;
208 
209  if (assHits.size()+hits.size()<maxfound) return;
210 
211  // loop over the remaining hits
212  for (hitCont::const_iterator hit=hits.begin(); hit!=hits.end(); ++hit) {
213 
214 // if (debug) {
215 // int nHits=assHits.size()+1;
216 // cout << " Trying B: " << **hit<< " wire: " << (*hit)->id() << endl;
217 // printPattern(assHits,*hit);
218 // }
219 
220  assHits.push_back(DTSegmentCand::AssPoint(*hit, DTEnums::Left));
221  std::unique_ptr<DTSegmentCand> left_seg = fitWithT0(sl,assHits, chi2l, t0_corrl,0);
222  assHits.pop_back();
223 // if (debug)
224 // cout << " Left: t0= " << t0_corrl << " chi2/nHits= " << chi2l << "/" << nHits << " ok: " << left_ok << endl;
225 
226  assHits.push_back(DTSegmentCand::AssPoint(*hit, DTEnums::Right));
227  std::unique_ptr<DTSegmentCand> right_seg = fitWithT0(sl,assHits, chi2r, t0_corrr,0);
228  assHits.pop_back();
229 // if (debug)
230 // cout << " Right: t0= " << t0_corrr << " chi2/nHits= " << chi2r << "/" << nHits << " ok: " << right_ok << endl;
231 
232  bool left_ok=(left_seg)?true:false;
233  bool right_ok=(right_seg)?true:false;
234 
235  if (!left_ok && !right_ok) continue;
236 
237  foundSomething = true;
238 
239  // prepare the hit set for the next search, start from the other side
240  hitCont hitsForFit;
241  for (hitCont::const_iterator tmpHit=hit+1; tmpHit!=hits.end(); tmpHit++)
242  if (geometryFilter((*tmpHit)->id(),(*hit)->id())) hitsForFit.push_back(*tmpHit);
243 
244  reverse(hitsForFit.begin(),hitsForFit.end());
245 
246  // choose only one - left or right
247  if (assHits.size()>3 && left_ok && right_ok) {
248  if (chi2l<chi2r-0.1) right_ok=false; else
249  if (chi2r<chi2l-0.1) left_ok=false;
250  }
251  if (left_ok) {
252  assHits.push_back(DTSegmentCand::AssPoint(*hit, DTEnums::Left));
253  addHits(sl,assHits,hitsForFit,result);
254  assHits.pop_back();
255  }
256  if (right_ok) {
257  assHits.push_back(DTSegmentCand::AssPoint(*hit, DTEnums::Right));
258  addHits(sl,assHits,hitsForFit,result);
259  assHits.pop_back();
260  }
261  }
262 
263  if (foundSomething) return;
264  // if we didn't find any new hits compatible with the current candidate, we proceed to check and store the candidate
265 
266  // If we already have a segment with more hits from this hit pair - don't save this one.
267  if (assHits.size()<maxfound) return;
268 
269  // Check if semgent Ok, calculate chi2
270  std::unique_ptr<DTSegmentCand> seg = fitWithT0(sl,assHits, chi2l, t0_corrl,debug);
271  if (!seg) return;
272 
273  if (!seg->good()) {
274 // if (debug) cout << " Segment not good() - skipping" << endl;
275  return;
276  }
277 
278  if (assHits.size()>maxfound) maxfound = assHits.size();
279  if (debug) cout << endl << " Seg t0= " << t0_corrl << endl << *seg<< endl;
280 
281  if (checkDoubleCandidates(result,seg.get())) {
282  result.push_back(seg.release());
283  if (debug) cout << " Result is now " << result.size() << endl;
284  } else {
285  if (debug) cout << " Exists - skipping" << endl;
286  }
287 }
288 
289 
290 bool
292 {
293 // return true;
294 
295  const int layerLowerCut[4]={0,-1,-2,-2};
296  const int layerUpperCut[4]={0, 2, 2, 3};
297 // const int layerLowerCut[4]={0,-2,-4,-5};
298 // const int layerUpperCut[4]={0, 3, 4, 6};
299 
300  // deal only with hits that are in the same SL
301  if (first.layerId().superlayerId().superLayer()!=second.layerId().superlayerId().superLayer())
302  return true;
303 
304  int deltaLayer=abs(first.layerId().layer()-second.layerId().layer());
305 
306  // drop hits in the same layer
307  if (!deltaLayer) return false;
308 
309  // protection against unexpected layer numbering
310  if (deltaLayer>3) {
311  cout << "*** WARNING! DT Layer numbers differ by more than 3! for hits: " << endl;
312  cout << " " << first << endl;
313  cout << " " << second << endl;
314  return false;
315  }
316 
317  // accept only hits in cells "not too far away"
318  int deltaWire=first.wire()-second.wire();
319  if (second.layerId().layer()%2==0) deltaWire=-deltaWire; // yet another trick to get it right...
320  if ((deltaWire<=layerLowerCut[deltaLayer]) || (deltaWire>=layerUpperCut[deltaLayer])) return false;
321 
322  return true;
323 }
324 
325 
326 std::unique_ptr<DTSegmentCand>
327 DTMeantimerPatternReco::fitWithT0(const DTSuperLayer* sl, const vector<DTSegmentCand::AssPoint> &assHits, double &chi2, double &t0_corr, const bool fitdebug)
328 {
329  // create a DTSegmentCand
330  DTSegmentCand::AssPointCont pointsSet;
331  pointsSet.insert(assHits.begin(),assHits.end());
332  std::unique_ptr<DTSegmentCand> seg(new DTSegmentCand(pointsSet,sl));
333 
334  // perform the 3 parameter fit on the segment candidate
335  theUpdator->fit(seg.get(),1,fitdebug);
336 
337  chi2 = seg->chi2();
338  t0_corr = seg->t0();
339 
340  // sanity check - drop segment candidates with a failed fit
341  // for a 3-par fit this includes segments with hits after the calculated t0 correction ending up
342  // beyond the chamber walls or on the other side of the wire
343  if (chi2==-1.) return nullptr;
344 
345  // at this point we keep all 3-hit segments that passed the above check
346  if (assHits.size()==3) return seg;
347 
348  // for segments with no t0 information we impose a looser chi2 cut
349  if (t0_corr==0) {
350  if (chi2<100.) return seg;
351  else return nullptr;
352  }
353 
354  // cut on chi2/ndof of the segment candidate
355  if ((chi2/(assHits.size()-3)<theMaxChi2)) return seg;
356  else return nullptr;
357 }
358 
359 
360 
361 bool
362 DTMeantimerPatternReco::checkDoubleCandidates(vector<DTSegmentCand*>& cands,
363  DTSegmentCand* seg) {
364  for (vector<DTSegmentCand*>::iterator cand=cands.begin();
365  cand!=cands.end(); ++cand) {
366  if (*(*cand)==*seg) return false;
367  if (((*cand)->nHits()>=seg->nHits()) && ((*cand)->chi2ndof()<seg->chi2ndof()))
368  if ((*cand)->nSharedHitPairs(*seg)>int(seg->nHits()-2)) return false;
369  }
370  return true;
371 }
372 
373 
374 
375 void
376 DTMeantimerPatternReco::printPattern( vector<DTSegmentCand::AssPoint>& assHits, const DTHitPairForFit* hit) {
377 
378  char mark[26]={". . . . . . . . . . . . "};
379  int wire[12]={0,0,0,0,0,0,0,0,0,0,0,0};
380 
381  for (vector<DTSegmentCand::AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
382  int lay = (((*assHit).first)->id().superlayerId().superLayer()-1)*4 + ((*assHit).first)->id().layerId().layer()-1;
383  wire[lay]= ((*assHit).first)->id().wire();
384  if ((*assHit).second==DTEnums::Left) mark[lay*2]='L'; else mark[lay*2]='R';
385  }
386 
387  int lay = ((*hit).id().superlayerId().superLayer()-1)*4 + (*hit).id().layerId().layer() - 1;
388  wire[lay]= (*hit).id().wire();
389  mark[lay*2]='*';
390 
391  cout << " " << mark << endl << " ";
392  for (int i=0; i<12; i++) if (wire[i]) cout << setw(2) << wire[i]; else cout << " ";
393  cout << endl;
394 }
std::pair< std::shared_ptr< DTHitPairForFit >, DTEnums::DTCellSide > AssPoint
Definition: DTSegmentCand.h:38
std::vector< std::shared_ptr< DTHitPairForFit > > initHits(const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual double chi2ndof() const
the chi2/NDOF of the fit
Definition: DTSegmentCand.h:67
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:52
DTCellSide
Which side of the DT cell.
Definition: DTEnums.h:15
virtual ~DTMeantimerPatternReco()
Destructor.
int layer() const
Return the layer number.
Definition: DTLayerId.h:53
std::vector< DTSegmentCand * > clean(const std::vector< DTSegmentCand * > &inputCands) const
do the cleaning
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:59
std::vector< std::shared_ptr< DTHitPairForFit > > hitCont
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:36
void push_back(D *&d)
Definition: OwnVector.h:280
std::vector< DTSegmentCand * > buildSegments(const DTSuperLayer *sl, const std::vector< std::shared_ptr< DTHitPairForFit >> &hits)
virtual void setES(const edm::EventSetup &setup)
void setES(const edm::EventSetup &setup)
set the setup
std::unique_ptr< DTSegmentCand > fitWithT0(const DTSuperLayer *sl, const std::vector< DTSegmentCand::AssPoint > &assHits, double &chi2, double &t0_corr, const bool fitdebug)
tuple result
Definition: query.py:137
tuple IP
Definition: listHistos.py:76
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int superLayer() const
Return the superlayer number.
bool first
Definition: L1TdeRCT.cc:75
int mark
Definition: cuy.py:1009
unsigned int id
int wire() const
Return the wire number.
Definition: DTWireId.h:56
void printPattern(std::vector< DTSegmentCand::AssPoint > &assHits, const DTHitPairForFit *hit)
std::set< AssPoint, AssPointLessZ > AssPointCont
Definition: DTSegmentCand.h:40
void addHits(const DTSuperLayer *sl, std::vector< DTSegmentCand::AssPoint > &assHits, const std::vector< std::shared_ptr< DTHitPairForFit >> &hits, std::vector< DTSegmentCand * > &result)
virtual unsigned int nHits() const
Definition: DTSegmentCand.h:61
hitCont::const_iterator hitIter
const T & get() const
Definition: EventSetup.h:55
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
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:62
DTMeantimerPatternReco(const edm::ParameterSet &pset)
Constructor.
tuple cout
Definition: gather_cfg.py:121
void update(DTRecSegment4D *seg, const bool calcT0, bool allow3par) const
recompute hits position and refit the segment4D
list pairs
sort tag files by run number
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
bool fit(DTSegmentCand *seg, bool allow3par, const bool fitdebug) const
edm::ESHandle< DTGeometry > theDTGeometry