CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCSegAlgoDF.cc
Go to the documentation of this file.
1 
7 #include "CSCSegAlgoDF.h"
8 
10 // For clhep Matrix::solve
12 
14 
17 
19 
21 #include "CSCSegAlgoShowering.h"
22 
23 #include <algorithm>
24 #include <cmath>
25 #include <iostream>
26 #include <string>
27 
28 
29 /* Constructor
30  *
31  */
32 CSCSegAlgoDF::CSCSegAlgoDF(const edm::ParameterSet& ps) : CSCSegmentAlgorithm(ps), myName("CSCSegAlgoDF") {
33 
34  debug = ps.getUntrackedParameter<bool>("CSCSegmentDebug");
35  minLayersApart = ps.getParameter<int>("minLayersApart");
36  minHitsPerSegment = ps.getParameter<int>("minHitsPerSegment");
37  dRPhiFineMax = ps.getParameter<double>("dRPhiFineMax");
38  dPhiFineMax = ps.getParameter<double>("dPhiFineMax");
39  tanThetaMax = ps.getParameter<double>("tanThetaMax");
40  tanPhiMax = ps.getParameter<double>("tanPhiMax");
41  chi2Max = ps.getParameter<double>("chi2Max");
42  preClustering = ps.getUntrackedParameter<bool>("preClustering");
43  minHitsForPreClustering= ps.getParameter<int>("minHitsForPreClustering");
44  nHitsPerClusterIsShower= ps.getParameter<int>("nHitsPerClusterIsShower");
45  Pruning = ps.getUntrackedParameter<bool>("Pruning");
46  maxRatioResidual = ps.getParameter<double>("maxRatioResidualPrune");
47 
49  showering_ = new CSCSegAlgoShowering( ps );
50 }
51 
52 
53 /* Destructor
54  *
55  */
57  delete preCluster_;
58  delete showering_;
59 }
60 
61 
62 /* run
63  *
64  */
65 std::vector<CSCSegment> CSCSegAlgoDF::run(const CSCChamber* aChamber, const ChamberHitContainer& rechits) {
66 
67  // Store chamber info in temp memory
68  theChamber = aChamber;
69 
70  int nHits = rechits.size();
71 
72  // Segments prior to pruning
73  std::vector<CSCSegment> segments_temp;
74 
75  if ( preClustering && nHits > minHitsForPreClustering ) {
76  // This is where the segment origin is in the chamber on avg.
77  std::vector<CSCSegment> testSegments;
78  std::vector<ChamberHitContainer> clusteredHits = preCluster_->clusterHits(theChamber, rechits);
79  // loop over the found clusters:
80  for (std::vector<ChamberHitContainer>::iterator subrechits = clusteredHits.begin(); subrechits != clusteredHits.end(); ++subrechits ) {
81  // build the subset of segments:
82  std::vector<CSCSegment> segs = buildSegments( (*subrechits) );
83  // add the found subset of segments to the collection of all segments in this chamber:
84  segments_temp.insert( segments_temp.end(), segs.begin(), segs.end() );
85  }
86  } else {
87  std::vector<CSCSegment> segs = buildSegments( rechits );
88  // add the found subset of segments to the collection of all segments in this chamber:
89  segments_temp.insert( segments_temp.end(), segs.begin(), segs.end() );
90  }
91 
92  return segments_temp;
93 }
94 
95 
96 /* This builds segments by first creating proto-segments from at least 3 hits.
97  * We intend to try all possible pairs of hits to start segment building. 'All possible'
98  * means each hit lies on different layers in the chamber. Once a hit has been assigned
99  * to a segment, we don't consider it again, THAT IS, FOR THE FIRST PASS ONLY !
100  * In fact, this is one of the possible flaw with the SK algorithms as it sometimes manages
101  * to build segments with the wrong starting points. In the DF algorithm, the endpoints
102  * are tested as the best starting points in a 2nd loop.
103  *
104  * Also, only a certain muonsPerChamberMax maximum number of segments can be produced in the chamber
105  */
106 std::vector<CSCSegment> CSCSegAlgoDF::buildSegments(const ChamberHitContainer& _rechits) {
107 
108  ChamberHitContainer rechits = _rechits;
109  // Clear buffer for segment vector
110  std::vector<CSCSegment> segmentInChamber;
111  segmentInChamber.clear();
112 
113  unsigned nHitInChamber = rechits.size();
114  if ( nHitInChamber < 3 ) return segmentInChamber;
115 
116  LayerIndex layerIndex( nHitInChamber );
117 
118  unsigned nLayers = 0;
119  int old_layer = -1;
120  for ( unsigned int i = 0; i < nHitInChamber; i++ ) {
121  int this_layer = rechits[i]->cscDetId().layer();
122  layerIndex[i] = this_layer;
123  if ( this_layer != old_layer ) {
124  old_layer = this_layer;
125  nLayers++;
126  }
127  }
128 
129  if ( nLayers < 3 ) return segmentInChamber;
130 
131  double z1 = theChamber->layer(1)->position().z();
132  double z6 = theChamber->layer(6)->position().z();
133 
134  if ( z1 > 0. ) {
135  if ( z1 > z6 ) {
136  reverse( layerIndex.begin(), layerIndex.end() );
137  reverse( rechits.begin(), rechits.end() );
138  }
139  }
140  else if ( z1 < 0. ) {
141  if ( z1 < z6 ) {
142  reverse( layerIndex.begin(), layerIndex.end() );
143  reverse( rechits.begin(), rechits.end() );
144  }
145  }
146 
147  // Showering muon
148  if ( preClustering && int(nHitInChamber) > nHitsPerClusterIsShower && nLayers > 2 ) {
149  CSCSegment segShower = showering_->showerSeg(theChamber, rechits);
150 
151  // Make sure have at least 3 hits...
152  if ( segShower.nRecHits() < 3 ) return segmentInChamber;
153 
154  segmentInChamber.push_back(segShower);
155 
156  return segmentInChamber;
157  }
158 
159 
160  // Initialize flags that a given hit has been allocated to a segment
161  BoolContainer used_ini(rechits.size(), false);
162  usedHits = used_ini;
163 
164  ChamberHitContainerCIt ib = rechits.begin();
165  ChamberHitContainerCIt ie = rechits.end();
166 
167  // Now Loop over hits within the chamber to find 1st seed for segment building
168  for ( ChamberHitContainerCIt i1 = ib; i1 < ie; ++i1 ) {
169  if ( usedHits[i1-ib] ) continue;
170 
171  const CSCRecHit2D* h1 = *i1;
172  int layer1 = layerIndex[i1-ib];
173  const CSCLayer* l1 = theChamber->layer(layer1);
174  GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
175  LocalPoint lp1 = theChamber->toLocal(gp1);
176 
177  // Loop over hits backward to find 2nd seed for segment building
178  for ( ChamberHitContainerCIt i2 = ie-1; i2 > ib; --i2 ) {
179 
180  if ( usedHits[i2-ib] ) continue; // Hit has been used already
181 
182  int layer2 = layerIndex[i2-ib];
183  if ( (layer2 - layer1) < minLayersApart ) continue;
184 
185  const CSCRecHit2D* h2 = *i2;
186  const CSCLayer* l2 = theChamber->layer(layer2);
187  GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
188  LocalPoint lp2 = theChamber->toLocal(gp2);
189 
190  // Clear proto segment so it can be (re)-filled
191  protoSegment.clear();
192 
193  // localPosition is position of hit wrt layer (so local z = 0)
195 
196  // We want hit wrt chamber (and local z will be != 0)
197  float dz = gp2.z()-gp1.z();
198  protoSlope_u = (lp2.x() - lp1.x())/dz ;
199  protoSlope_v = (lp2.y() - lp1.y())/dz ;
200 
201  // Test if entrance angle is roughly pointing towards IP
202  if (fabs(protoSlope_v) > tanThetaMax) continue;
203  if (fabs(protoSlope_u) > tanPhiMax ) continue;
204 
205  protoSegment.push_back(h1);
206  protoSegment.push_back(h2);
207 
208  // Try adding hits to proto segment
209  tryAddingHitsToSegment(rechits, i1, i2, layerIndex);
210 
211  // Check no. of hits on segment to see if segment is large enough
212  bool segok = true;
213  unsigned iadd = 0;
214 
215  if (protoSegment.size() < minHitsPerSegment+iadd) segok = false;
216 
217  if ( Pruning && segok ) pruneFromResidual();
218 
219  // Check if segment satisfies chi2 requirement
220  if (protoChi2 > chi2Max) segok = false;
221 
222  if ( segok ) {
223 
224  // Fill segment properties
225 
226  // Local direction
227  double dz = 1./sqrt(1. + protoSlope_u*protoSlope_u + protoSlope_v*protoSlope_v);
228  double dx = dz * protoSlope_u;
229  double dy = dz * protoSlope_v;
230  LocalVector localDir(dx,dy,dz);
231 
232  // localDir may need sign flip to ensure it points outward from IP
233  double globalZpos = ( theChamber->toGlobal( protoIntercept ) ).z();
234  double globalZdir = ( theChamber->toGlobal( localDir ) ).z();
235  double directionSign = globalZpos * globalZdir;
236  protoDirection = (directionSign * localDir).unit();
237 
238  // Error matrix
239  AlgebraicSymMatrix protoErrors = calculateError();
240  // but reorder components to match what's required by TrackingRecHit interface
241  // i.e. slopes first, then positions
242  flipErrors( protoErrors );
243 
245 
246  segmentInChamber.push_back(temp);
247 
248  if (nHitInChamber-protoSegment.size() < 3) return segmentInChamber;
249  if (segmentInChamber.size() > 4) return segmentInChamber;
250 
251  // Flag used hits
252  flagHitsAsUsed(rechits);
253  }
254  }
255  }
256  return segmentInChamber;
257 }
258 
259 
260 /* Method tryAddingHitsToSegment
261  *
262  * Look at left over hits and try to add them to proto segment by looking how far they
263  * are from the segment in terms of the hit error matrix (so how many sigmas away).
264  *
265  */
267  const ChamberHitContainerCIt i1,
268  const ChamberHitContainerCIt i2,
269  const LayerIndex& layerIndex ) {
270 
271 /* Iterate over the layers with hits in the chamber
272  * Skip the layers containing the segment endpoints on first pass, but then
273  * try hits on layer containing the segment starting points on 2nd pass
274  * if segment has >2 hits. Once a hit is added to a layer, don't replace it
275  * until 2nd iteration
276  */
277 
278  ChamberHitContainerCIt ib = rechits.begin();
279  ChamberHitContainerCIt ie = rechits.end();
280  closeHits.clear();
281 
282  for ( ChamberHitContainerCIt i = ib; i != ie; ++i ) {
283 
284  if (i == i1 || i == i2 ) continue;
285  if ( usedHits[i-ib] ) continue; // Don't use hits already part of a segment.
286 
287  const CSCRecHit2D* h = *i;
288  int layer = layerIndex[i-ib];
289  int layer1 = layerIndex[i1-ib];
290  int layer2 = layerIndex[i2-ib];
291 
292  // Low multiplicity case
293  if (rechits.size() < 9) {
294  if ( isHitNearSegment( h ) ) {
295  if ( !hasHitOnLayer(layer) ) {
296  addHit(h, layer);
297  } else {
298  closeHits.push_back(h);
299  }
300  }
301 
302  // High multiplicity case
303  } else {
304  if ( isHitNearSegment( h ) ) {
305  if ( !hasHitOnLayer(layer) ) {
306  addHit(h, layer);
308  // Don't change the starting points at this stage !!!
309  } else {
310  closeHits.push_back(h);
311  if (layer != layer1 && layer != layer2 ) compareProtoSegment(h, layer);
312  }
313  }
314  }
315  }
316 
317  if ( int(protoSegment.size()) < 3) return;
318 
320 
321  // 2nd pass to remove biases
322  // This time, also consider changing the endpoints
323  for ( ChamberHitContainerCIt i = closeHits.begin() ; i != closeHits.end(); ++i ) {
324  const CSCRecHit2D* h = *i;
325  int layer = (*i)->cscDetId().layer();
326  compareProtoSegment(h, layer);
327  }
328 
329 }
330 
331 
332 /* isHitNearSegment
333  *
334  * Compare rechit with expected position from proto_segment
335  */
337 
338  const CSCLayer* layer = theChamber->layer(hit->cscDetId().layer());
339 
340  // hit phi position in global coordinates
341  GlobalPoint Hgp = layer->toGlobal(hit->localPosition());
342  double Hphi = Hgp.phi();
343  if (Hphi < 0.) Hphi += 2.*M_PI;
344  LocalPoint Hlp = theChamber->toLocal(Hgp);
345  double z = Hlp.z();
346 
347  double LocalX = protoIntercept.x() + protoSlope_u * z;
348  double LocalY = protoIntercept.y() + protoSlope_v * z;
349  LocalPoint Slp(LocalX, LocalY, z);
350  GlobalPoint Sgp = theChamber->toGlobal(Slp);
351  double Sphi = Sgp.phi();
352  if (Sphi < 0.) Sphi += 2.*M_PI;
353  double R = sqrt(Sgp.x()*Sgp.x() + Sgp.y()*Sgp.y());
354 
355  double deltaPhi = Sphi - Hphi;
356  if (deltaPhi > 2.*M_PI) deltaPhi -= 2.*M_PI;
357  if (deltaPhi < -2.*M_PI) deltaPhi += 2.*M_PI;
358  if (deltaPhi < 0.) deltaPhi = -deltaPhi;
359 
360  double RdeltaPhi = R * deltaPhi;
361 
362  if (RdeltaPhi < dRPhiFineMax && deltaPhi < dPhiFineMax ) return true;
363 
364  return false;
365 }
366 
367 
368 /* Method addHit
369  *
370  * Test if can add hit to proto segment. If so, try to add it.
371  *
372  */
373 bool CSCSegAlgoDF::addHit(const CSCRecHit2D* aHit, int layer) {
374 
375  // Return true if hit was added successfully and then parameters are updated.
376  // Return false if there is already a hit on the same layer, or insert failed.
377 
378  bool ok = true;
379 
380  // Test that we are not trying to add the same hit again
381  for ( ChamberHitContainer::const_iterator it = protoSegment.begin(); it != protoSegment.end(); it++ )
382  if ( aHit == (*it) ) return false;
383 
384  protoSegment.push_back(aHit);
385 
386  return ok;
387 }
388 
389 
390 /* Method updateParameters
391  *
392  * Perform a simple Least Square Fit on proto segment to determine slope and intercept
393  *
394  */
396 
397  // Compute slope from Least Square Fit
398  CLHEP::HepMatrix M(4,4,0);
399  CLHEP::HepVector B(4,0);
400 
401  ChamberHitContainer::const_iterator ih;
402 
403  for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
404 
405  const CSCRecHit2D& hit = (**ih);
406  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
407  GlobalPoint gp = layer->toGlobal(hit.localPosition());
408  LocalPoint lp = theChamber->toLocal(gp);
409 
410  double u = lp.x();
411  double v = lp.y();
412  double z = lp.z();
413 
414  // ptc: Covariance matrix of local errors
415  CLHEP::HepMatrix IC(2,2);
416  IC(1,1) = hit.localPositionError().xx();
417  IC(1,2) = hit.localPositionError().xy();
418  IC(2,2) = hit.localPositionError().yy();
419  IC(2,1) = IC(1,2); // since Cov is symmetric
420 
421  // ptc: Invert covariance matrix (and trap if it fails!)
422  int ierr = 0;
423  IC.invert(ierr); // inverts in place
424  if (ierr != 0) {
425  LogDebug("CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC << "\n";
426  }
427 
428  M(1,1) += IC(1,1);
429  M(1,2) += IC(1,2);
430  M(1,3) += IC(1,1) * z;
431  M(1,4) += IC(1,2) * z;
432  B(1) += u * IC(1,1) + v * IC(1,2);
433 
434  M(2,1) += IC(2,1);
435  M(2,2) += IC(2,2);
436  M(2,3) += IC(2,1) * z;
437  M(2,4) += IC(2,2) * z;
438  B(2) += u * IC(2,1) + v * IC(2,2);
439 
440  M(3,1) += IC(1,1) * z;
441  M(3,2) += IC(1,2) * z;
442  M(3,3) += IC(1,1) * z * z;
443  M(3,4) += IC(1,2) * z * z;
444  B(3) += ( u * IC(1,1) + v * IC(1,2) ) * z;
445 
446  M(4,1) += IC(2,1) * z;
447  M(4,2) += IC(2,2) * z;
448  M(4,3) += IC(2,1) * z * z;
449  M(4,4) += IC(2,2) * z * z;
450  B(4) += ( u * IC(2,1) + v * IC(2,2) ) * z;
451  }
452 
453  CLHEP::HepVector p = solve(M, B);
454 
455  // Update member variables
456  // Note that origin has local z = 0
457 
458  protoIntercept = LocalPoint(p(1), p(2), 0.);
459  protoSlope_u = p(3);
460  protoSlope_v = p(4);
461 
462 
463  // Determine Chi^2 for the proto wire segment
464 
465  double chsq = 0.;
466 
467  for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
468 
469  const CSCRecHit2D& hit = (**ih);
470  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
471  GlobalPoint gp = layer->toGlobal(hit.localPosition());
472  LocalPoint lp = theChamber->toLocal(gp);
473 
474  double u = lp.x();
475  double v = lp.y();
476  double z = lp.z();
477 
478  double du = protoIntercept.x() + protoSlope_u * z - u;
479  double dv = protoIntercept.y() + protoSlope_v * z - v;
480 
481  CLHEP::HepMatrix IC(2,2);
482  IC(1,1) = hit.localPositionError().xx();
483  IC(1,2) = hit.localPositionError().xy();
484  IC(2,2) = hit.localPositionError().yy();
485  IC(2,1) = IC(1,2);
486 
487  // Invert covariance matrix
488  int ierr = 0;
489  IC.invert(ierr);
490  if (ierr != 0) {
491  LogDebug("CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC << "\n";
492  }
493  chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
494  }
495  protoChi2 = chsq;
496 }
497 
498 
499 /* hasHitOnLayer
500  *
501  * Just make sure hit to be added to layer comes from different layer than those in proto segment
502  */
503 bool CSCSegAlgoDF::hasHitOnLayer(int layer) const {
504 
505  // Is there already a hit on this layer?
506  for ( ChamberHitContainerCIt it = protoSegment.begin(); it != protoSegment.end(); it++ )
507  if ( (*it)->cscDetId().layer() == layer ) return true;
508 
509  return false;
510 }
511 
512 
513 /* Method compareProtoSegment
514  *
515  * For hit coming from the same layer of an existing hit within the proto segment
516  * test if achieve better chi^2 by using this hit than the other
517  *
518  */
520 
521  // Store old segment first
522  double old_protoChi2 = protoChi2;
523  LocalPoint old_protoIntercept = protoIntercept;
524  float old_protoSlope_u = protoSlope_u;
525  float old_protoSlope_v = protoSlope_v;
526  LocalVector old_protoDirection = protoDirection;
527  ChamberHitContainer old_protoSegment = protoSegment;
528 
529 
530  // Try adding the hit to existing segment, and remove old one existing in same layer
531  ChamberHitContainer::iterator it;
532  for ( it = protoSegment.begin(); it != protoSegment.end(); ) {
533  if ( (*it)->cscDetId().layer() == layer ) {
534  it = protoSegment.erase(it);
535  } else {
536  ++it;
537  }
538  }
539  bool ok = addHit(h, layer);
540 
541  if (ok) updateParameters();
542 
543  if ( (protoChi2 > old_protoChi2) || ( !ok ) ) {
544  protoChi2 = old_protoChi2;
545  protoIntercept = old_protoIntercept;
546  protoSlope_u = old_protoSlope_u;
547  protoSlope_v = old_protoSlope_v;
548  protoDirection = old_protoDirection;
549  protoSegment = old_protoSegment;
550  }
551 }
552 
553 
554 /* Method flagHitsAsUsed
555  *
556  * Flag hits which have entered segment building so we don't reuse them.
557  * Also flag does which were very close to segment to reduce combinatorics
558  */
559 void CSCSegAlgoDF::flagHitsAsUsed(const ChamberHitContainer& rechitsInChamber) {
560 
561  // Flag hits on segment as used
562  ChamberHitContainerCIt ib = rechitsInChamber.begin();
563  ChamberHitContainerCIt hi, iu;
564 
565  for ( hi = protoSegment.begin(); hi != protoSegment.end(); ++hi ) {
566  for ( iu = ib; iu != rechitsInChamber.end(); ++iu ) {
567  if (*hi == *iu) usedHits[iu-ib] = true;
568  }
569  }
570  // Don't reject hits marked as "nearby" for now.
571  // So this is bypassed at all times for now !!!
572  // Perhaps add back to speed up algorithm some more
573  if (closeHits.size() > 0) return;
574  for ( hi = closeHits.begin(); hi != closeHits.end(); ++hi ) {
575  for ( iu = ib; iu != rechitsInChamber.end(); ++iu ) {
576  if (*hi == *iu) usedHits[iu-ib] = true;
577  }
578  }
579 
580 }
581 
582 
584 
585  std::vector<const CSCRecHit2D*>::const_iterator it;
586  int nhits = protoSegment.size();
587  AlgebraicSymMatrix matrix(2*nhits, 0);
588  int row = 0;
589 
590  for (it = protoSegment.begin(); it != protoSegment.end(); ++it) {
591 
592  const CSCRecHit2D& hit = (**it);
593  ++row;
594  matrix(row, row) = hit.localPositionError().xx();
595  matrix(row, row+1) = hit.localPositionError().xy();
596  ++row;
597  matrix(row, row-1) = hit.localPositionError().xy();
598  matrix(row, row) = hit.localPositionError().yy();
599  }
600  int ierr;
601  matrix.invert(ierr);
602  return matrix;
603 }
604 
605 
606 CLHEP::HepMatrix CSCSegAlgoDF::derivativeMatrix() const {
607 
608  ChamberHitContainer::const_iterator it;
609  int nhits = protoSegment.size();
610  CLHEP::HepMatrix matrix(2*nhits, 4);
611  int row = 0;
612 
613  for(it = protoSegment.begin(); it != protoSegment.end(); ++it) {
614 
615  const CSCRecHit2D& hit = (**it);
616  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
617  GlobalPoint gp = layer->toGlobal(hit.localPosition());
618  LocalPoint lp = theChamber->toLocal(gp);
619  float z = lp.z();
620  ++row;
621  matrix(row, 1) = 1.;
622  matrix(row, 3) = z;
623  ++row;
624  matrix(row, 2) = 1.;
625  matrix(row, 4) = z;
626  }
627  return matrix;
628 }
629 
630 /* calculateError
631  *
632  */
634 
637 
638  // (AT W A)^-1
639  // from http://www.phys.ufl.edu/~avery/fitting.html, part I
640  int ierr;
641  AlgebraicSymMatrix result = weights.similarityT(A);
642  result.invert(ierr);
643 
644  // blithely assuming the inverting never fails...
645  return result;
646 }
647 
649 
650  // The CSCSegment needs the error matrix re-arranged
651 
652  AlgebraicSymMatrix hold( a );
653 
654  // errors on slopes into upper left
655  a(1,1) = hold(3,3);
656  a(1,2) = hold(3,4);
657  a(2,1) = hold(4,3);
658  a(2,2) = hold(4,4);
659 
660  // errors on positions into lower right
661  a(3,3) = hold(1,1);
662  a(3,4) = hold(1,2);
663  a(4,3) = hold(2,1);
664  a(4,4) = hold(2,2);
665 
666  // off-diagonal elements remain unchanged
667 
668 }
669 
670 
671 // Try to clean up segments by quickly looking at residuals
673 
674  // Only prune if have at least 5 hits
675  if ( protoSegment.size() < 5 ) return ;
676 
677 
678  // Now Study residuals
679 
680  float maxResidual = 0.;
681  float sumResidual = 0.;
682  int nHits = 0;
683  int badIndex = -1;
684  int j = 0;
685 
686 
687  ChamberHitContainer::const_iterator ih;
688 
689  for ( ih = protoSegment.begin(); ih != protoSegment.end(); ++ih ) {
690  const CSCRecHit2D& hit = (**ih);
691  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
692  GlobalPoint gp = layer->toGlobal(hit.localPosition());
693  LocalPoint lp = theChamber->toLocal(gp);
694 
695  double u = lp.x();
696  double v = lp.y();
697  double z = lp.z();
698 
699  double du = protoIntercept.x() + protoSlope_u * z - u;
700  double dv = protoIntercept.y() + protoSlope_v * z - v;
701 
702  float residual = sqrt(du*du + dv*dv);
703 
704  sumResidual += residual;
705  nHits++;
706  if ( residual > maxResidual ) {
707  maxResidual = residual;
708  badIndex = j;
709  }
710  j++;
711  }
712 
713  float corrAvgResidual = (sumResidual - maxResidual)/(nHits -1);
714 
715  // Keep all hits
716  if ( maxResidual/corrAvgResidual < maxRatioResidual ) return;
717 
718 
719  // Drop worse hit and recompute segment properties + fill
720 
721  ChamberHitContainer newProtoSegment;
722 
723  j = 0;
724  for ( ih = protoSegment.begin(); ih != protoSegment.end(); ++ih ) {
725  if ( j != badIndex ) newProtoSegment.push_back(*ih);
726  j++;
727  }
728 
729  protoSegment.clear();
730 
731  for ( ih = newProtoSegment.begin(); ih != newProtoSegment.end(); ++ih ) {
732  protoSegment.push_back(*ih);
733  }
734 
735  // Update segment parameters
737 
738 }
739 
740 
741 /*
742  * Order the hits such that 2nd one is closest in x,y to first seed hit in global coordinates
743  */
745  const ChamberHitContainerCIt i1,
746  const ChamberHitContainerCIt i2,
748  const LayerIndex& layerIndex ) {
749 
750  secondSeedHits.clear();
751 
752  //ChamberHitContainerCIt ib = rechits.begin();
753  ChamberHitContainerCIt ie = rechits.end();
754 
755  // int layer1 = layerIndex[i1-ib];
756  // int layer2 = layerIndex[i2-ib];
757 
758 
759  // Now fill vector of rechits closest to center of mass:
760  // secondSeedHitsIdx.clear() = 0;
761 
762  // Loop over all hits and find hit closest to 1st seed.
763  for ( ChamberHitContainerCIt i2 = ie-1; i2 > i1; --i2 ) {
764 
765 
766  }
767 
768 
769 }
#define LogDebug(id)
void tryAddingHitsToSegment(const ChamberHitContainer &rechitsInChamber, const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2, const LayerIndex &layerIndex)
Utility functions.
void flipErrors(AlgebraicSymMatrix &) const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
float xx() const
Definition: LocalError.h:24
ChamberHitContainer protoSegment
Definition: CSCSegAlgoDF.h:135
std::vector< const CSCRecHit2D * > ChamberHitContainer
Definition: CSCSegAlgoDF.h:55
std::vector< const CSCRecHit2D * >::const_iterator ChamberHitContainerCIt
Definition: CSCSegAlgoDF.h:56
void orderSecondSeed(GlobalPoint gp1, const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2, const ChamberHitContainer &rechits, const LayerIndex &layerIndex)
CSCDetId cscDetId() const
Definition: CSCRecHit2D.h:52
int ib
Definition: cuy.py:660
AlgebraicSymMatrix weightMatrix(void) const
CLHEP::HepMatrix derivativeMatrix(void) const
virtual ~CSCSegAlgoDF()
Destructor.
Definition: CSCSegAlgoDF.cc:56
CSCSegAlgoPreClustering * preCluster_
Definition: CSCSegAlgoDF.h:161
double_binary B
Definition: DDStreamer.cc:234
CSCSegment showerSeg(const CSCChamber *aChamber, const ChamberHitContainer &rechits)
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:52
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
ChamberHitContainer closeHits
Definition: CSCSegAlgoDF.h:133
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:67
bool addHit(const CSCRecHit2D *hit, int layer)
bool isHitNearSegment(const CSCRecHit2D *h) const
CSCSegAlgoDF(const edm::ParameterSet &ps)
Constructor.
Definition: CSCSegAlgoDF.cc:32
CSCSegAlgoShowering * showering_
Definition: CSCSegAlgoDF.h:162
int layer() const
Definition: CSCDetId.h:74
ChamberHitContainer secondSeedHits
Definition: CSCSegAlgoDF.h:136
float float float z
bool preClustering
Definition: CSCSegAlgoDF.h:145
std::vector< CSCSegment > buildSegments(const ChamberHitContainer &rechits)
std::vector< CSCSegment > run(const CSCChamber *aChamber, const ChamberHitContainer &rechits)
Definition: CSCSegAlgoDF.cc:65
float xy() const
Definition: LocalError.h:25
void updateParameters(void)
std::deque< bool > BoolContainer
Definition: CSCSegAlgoDF.h:57
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:46
int nRecHits() const
Definition: CSCSegment.h:67
CLHEP::HepMatrix AlgebraicMatrix
LocalError localPositionError() const
Definition: CSCRecHit2D.h:51
float yy() const
Definition: LocalError.h:26
LocalPoint protoIntercept
Definition: CSCSegAlgoDF.h:139
string unit
Definition: csvLumiCalc.py:46
T sqrt(T t)
Definition: SSEVec.h:48
void compareProtoSegment(const CSCRecHit2D *h, int layer)
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:39
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
std::vector< int > LayerIndex
Typedefs.
Definition: CSCSegAlgoDF.h:54
float maxRatioResidual
Definition: CSCSegAlgoDF.h:159
#define M_PI
int minHitsForPreClustering
Definition: CSCSegAlgoDF.h:146
float protoSlope_v
Definition: CSCSegAlgoDF.h:138
void flagHitsAsUsed(const ChamberHitContainer &rechitsInChamber)
float tanThetaMax
Definition: CSCSegAlgoDF.h:157
double dPhiFineMax
Definition: CSCSegAlgoDF.h:155
float protoSlope_u
Definition: CSCSegAlgoDF.h:137
std::vector< std::vector< const CSCRecHit2D * > > clusterHits(const CSCChamber *aChamber, const ChamberHitContainer &rechits)
clusterize
double dRPhiFineMax
Definition: CSCSegAlgoDF.h:154
double a
Definition: hdecay.h:121
void pruneFromResidual()
CLHEP::HepSymMatrix AlgebraicSymMatrix
double protoChi2
Definition: CSCSegAlgoDF.h:140
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
AlgebraicSymMatrix calculateError(void) const
const CSCChamber * theChamber
Definition: CSCSegAlgoDF.h:130
LocalVector protoDirection
Definition: CSCSegAlgoDF.h:141
LocalPoint localPosition() const
Definition: CSCRecHit2D.h:50
int minHitsPerSegment
Definition: CSCSegAlgoDF.h:152
int nHitsPerClusterIsShower
Definition: CSCSegAlgoDF.h:150
T x() const
Definition: PV3DBase.h:62
BoolContainer usedHits
Definition: CSCSegAlgoDF.h:131
bool hasHitOnLayer(int layer) const