CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCSegAlgoTC.cc
Go to the documentation of this file.
1 
8 #include "CSCSegAlgoTC.h"
9 
12 
13 // For clhep Matrix::solve
16 
19 
21 
22 #include <algorithm>
23 #include <cmath>
24 #include <iostream>
25 #include <string>
26 
28  myName("CSCSegAlgoTC") {
29 
30  debugInfo = ps.getUntrackedParameter<bool>("verboseInfo");
31 
32  dRPhiMax = ps.getParameter<double>("dRPhiMax");
33  dPhiMax = ps.getParameter<double>("dPhiMax");
34  dRPhiFineMax = ps.getParameter<double>("dRPhiFineMax");
35  dPhiFineMax = ps.getParameter<double>("dPhiFineMax");
36  chi2Max = ps.getParameter<double>("chi2Max");
37  chi2ndfProbMin = ps.getParameter<double>("chi2ndfProbMin");
38  minLayersApart = ps.getParameter<int>("minLayersApart");
39  SegmentSorting = ps.getParameter<int>("SegmentSorting");
40 
41  LogDebug("CSC") << myName << " has algorithm cuts set to: \n"
42  << "--------------------------------------------------------------------\n"
43  << "dRPhiMax = " << dRPhiMax << '\n'
44  << "dPhiMax = " << dPhiMax << '\n'
45  << "dRPhiFineMax = " << dRPhiFineMax << '\n'
46  << "dPhiFineMax = " << dPhiFineMax << '\n'
47  << "chi2Max = " << chi2Max << '\n'
48  << "chi2ndfProbMin = " << chi2ndfProbMin << '\n'
49  << "minLayersApart = " << minLayersApart << '\n'
50  << "SegmentSorting = " << SegmentSorting << std::endl;
51 }
52 
53 std::vector<CSCSegment> CSCSegAlgoTC::run(const CSCChamber* aChamber, const ChamberHitContainer& rechits) {
54  theChamber = aChamber;
55  return buildSegments(rechits);
56 }
57 
58 std::vector<CSCSegment> CSCSegAlgoTC::buildSegments(const ChamberHitContainer& _rechits) {
59 
60  // Reimplementation of original algorithm of CSCSegmentizer, Mar-06
61 
62  LogDebug("CSC") << "*********************************************";
63  LogDebug("CSC") << "Start segment building in the new chamber: " << theChamber->specs()->chamberTypeName();
64  LogDebug("CSC") << "*********************************************";
65  ChamberHitContainer rechits = _rechits;
66  LayerIndex layerIndex(rechits.size());
67 
68  for(unsigned int i = 0; i < rechits.size(); i++) {
69 
70  layerIndex[i] = rechits[i]->cscDetId().layer();
71  }
72 
73  double z1 = theChamber->layer(1)->position().z();
74  double z6 = theChamber->layer(6)->position().z();
75 
76  if ( z1 > 0. ) {
77  if ( z1 > z6 ) {
78  reverse(layerIndex.begin(), layerIndex.end());
79  reverse(rechits.begin(), rechits.end());
80  }
81  }
82  else if ( z1 < 0. ) {
83  if ( z1 < z6 ) {
84  reverse(layerIndex.begin(), layerIndex.end());
85  reverse(rechits.begin(), rechits.end());
86  }
87  }
88 
89  if (debugInfo) {
90  // dump after sorting
91  dumpHits(rechits);
92  }
93 
94  if (rechits.size() < 2) {
95  LogDebug("CSC") << myName << ": " << rechits.size() <<
96  " hit(s) in chamber is not enough to build a segment.\n";
97  return std::vector<CSCSegment>();
98  }
99 
100  // We have at least 2 hits. We intend to try all possible pairs of hits to start
101  // segment building. 'All possible' means each hit lies on different layers in the chamber.
102  // For now we don't care whether a hit has already been allocated to another segment;
103  // we'll sort that out after building all possible segments.
104 
105  // Choose first hit (as close to IP as possible) h1 and
106  // second hit (as far from IP as possible) h2
107  // To do this we iterate over hits in the chamber by layer - pick two layers.
108  // @@ Require the two layers are at least 3 layers apart. May need tuning?
109  // Then we iterate over hits within each of these layers and pick h1 and h2 from these.
110  // If they are 'close enough' we build an empty segment.
111  // Then try adding hits to this segment.
112 
113  // Define buffer for segments we build (at the end we'll sort them somehow, and remove
114  // those which share hits with better-quality segments.
115 
116 
117  std::vector<CSCSegment> segments;
118 
119  ChamberHitContainerCIt ib = rechits.begin();
120  ChamberHitContainerCIt ie = rechits.end();
121 
122  for (ChamberHitContainerCIt i1 = ib; i1 != ie; ++i1) {
123 
124  int layer1 = layerIndex[i1-ib];
125  const CSCRecHit2D* h1 = *i1;
126 
127  for (ChamberHitContainerCIt i2 = ie-1; i2 != i1; --i2) {
128 
129  int layer2 = layerIndex[i2-ib];
130 
131  if (abs(layer2 - layer1) < minLayersApart)
132  break;
133 
134  const CSCRecHit2D* h2 = *i2;
135 
136  if (areHitsCloseInLocalX(h1, h2) && areHitsCloseInGlobalPhi(h1, h2)) {
137 
138  proto_segment.clear();
139 
140  const CSCLayer* l1 = theChamber->layer(h1->cscDetId().layer());
141  GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
142  const CSCLayer* l2 = theChamber->layer(h2->cscDetId().layer());
143  GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
144  LogDebug("CSC") << "start new segment from hits " << "h1: " << gp1 << " - h2: " << gp2 << "\n";
145 
146  if (!addHit(h1, layer1)) {
147  LogDebug("CSC") << " failed to add hit h1\n";
148  continue;
149  }
150 
151  if (!addHit(h2, layer2)) {
152  LogDebug("CSC") << " failed to add hit h2\n";
153  continue;
154  }
155 
156  tryAddingHitsToSegment(rechits, i1, i2); // changed seg
157 
158  // if a segment has been found push back it into the segment collection
159  if (proto_segment.empty()) {
160 
161  LogDebug("CSC") << "No segment has been found !!!\n";
162  }
163  else {
164 
165  // calculate error matrix
166  AlgebraicSymMatrix error_matrix = calculateError();
167 
168  // but reorder components to match what's required by TrackingRecHit interface
169  // i.e. slopes first, then positions
170  flipErrors( error_matrix );
171 
172  candidates.push_back(proto_segment);
173  origins.push_back(theOrigin);
174  directions.push_back(theDirection);
175  errors.push_back(error_matrix);
176  chi2s.push_back(theChi2);
177  LogDebug("CSC") << "Found a segment !!!\n";
178  }
179  }
180  }
181  }
182 
183  // We've built all possible segments. Now pick the best, non-overlapping subset.
184  pruneTheSegments(rechits);
185 
186  // Copy the selected proto segments into the CSCSegment vector
187  for(unsigned int i=0; i < candidates.size(); i++) {
188 
190  segments.push_back(temp);
191  }
192 
193  candidates.clear();
194  origins.clear();
195  directions.clear();
196  errors.clear();
197  chi2s.clear();
198 
199  // Give the segments to the CSCChamber
200  return segments;
201 }
202 
204  const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2) {
205 
206  // Iterate over the layers with hits in the chamber
207  // Skip the layers containing the segment endpoints
208  // Test each hit on the other layers to see if it is near the segment
209  // If it is, see whether there is already a hit on the segment from the same layer
210  // - if so, and there are more than 2 hits on the segment, copy the segment,
211  // replace the old hit with the new hit. If the new segment chi2 is better
212  // then replace the original segment with the new one (by swap)
213  // - if not, copy the segment, add the hit. If the new chi2/dof is still satisfactory
214  // then replace the original segment with the new one (by swap)
215 
216  ChamberHitContainerCIt ib = rechits.begin();
217  ChamberHitContainerCIt ie = rechits.end();
218 
219  for (ChamberHitContainerCIt i = ib; i != ie; ++i) {
220 
221  if ( i == i1 || i == i2 )
222  continue;
223 
224  int layer = (*i)->cscDetId().layer();
225  const CSCRecHit2D* h = *i;
226 
227  if (isHitNearSegment(h)) {
228 
229  const CSCLayer* l1 = theChamber->layer(h->cscDetId().layer());
230  GlobalPoint gp1 = l1->toGlobal(h->localPosition());
231  LogDebug("CSC") << " hit at global " << gp1 << " is near segment\n.";
232 
233  if (hasHitOnLayer(layer)) {
234  if (proto_segment.size() <= 2) {
235  LogDebug("CSC") << " " << proto_segment.size()
236  << " hits on segment...skip hit on same layer.\n";
237  continue;
238  }
239 
240  compareProtoSegment(h, layer);
241  }
242  else
243  increaseProtoSegment(h, layer);
244  } // h & seg close
245  } // i
246 }
247 
248 bool CSCSegAlgoTC::addHit(const CSCRecHit2D* aHit, int layer) {
249 
250  // Return true if hit was added successfully
251  // (and then parameters are updated).
252  // Return false if there is already a hit on the same layer, or insert failed.
253  bool ok = true;
254  ChamberHitContainer::const_iterator it;
255 
256  for(it = proto_segment.begin(); it != proto_segment.end(); it++)
257  if (((*it)->cscDetId().layer() == layer) && (aHit != *it))
258  return false;
259 
260  if (ok) {
261  proto_segment.push_back(aHit);
263  }
264  return ok;
265 }
266 
267 bool CSCSegAlgoTC::replaceHit(const CSCRecHit2D* h, int layer) {
268 
269  // replace a hit from a layer
270  ChamberHitContainer::iterator it;
271  for (it = proto_segment.begin(); it != proto_segment.end();) {
272  if ((*it)->cscDetId().layer() == layer)
273  it = proto_segment.erase(it);
274  else
275  ++it;
276  }
277 
278  return addHit(h, layer);
279 }
280 
282 
283  // Note that we need local position of a RecHit w.r.t. the CHAMBER
284  // and the RecHit itself only knows its local position w.r.t.
285  // the LAYER, so need to explicitly transform to global position.
286 
287  // no. of hits in the RecHitsOnSegment
288  // By construction this is the no. of layers with hitsna parte รจ da inserirsi tra le Contrade aperte ad accettare quello che
289 
290 
291  // since we allow just one hit per layer in a segment.
292 
293  int nh = proto_segment.size();
294 
295  // First hit added to a segment must always fail here
296  if (nh < 2)
297  return;
298 
299  if (nh == 2) {
300 
301  // Once we have two hits we can calculate a straight line
302  // (or rather, a straight line for each projection in xz and yz.)
303  ChamberHitContainer::const_iterator ih = proto_segment.begin();
304  int il1 = (*ih)->cscDetId().layer();
305  const CSCRecHit2D& h1 = (**ih);
306  ++ih;
307  int il2 = (*ih)->cscDetId().layer();
308  const CSCRecHit2D& h2 = (**ih);
309 
310  //@@ Skip if on same layer, but should be impossible
311  if (il1 == il2)
312  return;
313 
314  const CSCLayer* layer1 = theChamber->layer(il1);
315  const CSCLayer* layer2 = theChamber->layer(il2);
316 
317  GlobalPoint h1glopos = layer1->toGlobal(h1.localPosition());
318  GlobalPoint h2glopos = layer2->toGlobal(h2.localPosition());
319 
320  // localPosition is position of hit wrt layer (so local z = 0)
321  theOrigin = h1.localPosition();
322 
323  // We want hit wrt chamber (and local z will be != 0)
324  LocalPoint h1pos = theChamber->toLocal(h1glopos); // FIX !!
325  LocalPoint h2pos = theChamber->toLocal(h2glopos); // FIX !!
326 
327  float dz = h2pos.z()-h1pos.z();
328  uz = (h2pos.x()-h1pos.x())/dz ;
329  vz = (h2pos.y()-h1pos.y())/dz ;
330 
331  theChi2 = 0.;
332  }
333  else if (nh > 2) {
334 
335  // When we have more than two hits then we can fit projections to straight lines
336  fitSlopes();
337  fillChiSquared();
338  } // end of 'if' testing no. of hits
339 
341 }
342 
344 
345  // Update parameters of fit
346  // ptc 13-Aug-02: This does a linear least-squares fit
347  // to the hits associated with the segment, in the z projection.
348 
349  // In principle perhaps one could fit the strip and wire
350  // measurements (u, v respectively), to
351  // u = u0 + uz * z
352  // v = v0 + vz * z
353  // where u0, uz, v0, vz are the parameters resulting from the fit.
354  // But what is actually done is fit to the local x, y coordinates
355  // of the RecHits. However the strip measurement controls the precision
356  // of x, and the wire measurement controls that of y.
357  // Precision in local coordinate:
358  // u (strip, sigma~200um), v (wire, sigma~1cm)
359 
360  // I have verified that this code agrees with the formulation given
361  // on p246-247 of 'Data analysis techniques for high-energy physics
362  // experiments' by Bock, Grote, Notz & Regler, and that on p111-113
363  // of 'Statistics' by Barlow.
364 
365  // Formulate the matrix equation representing the least-squares fit
366  // We have a vector of measurements m, which is a 2n x 1 dim matrix
367  // The transpose mT is (u1, v1, u2, v2, ..., un, vn)
368  // where ui is the strip-associated measurement and vi is the
369  // wire-associated measurement for a given RecHit i.
370  // The fit is to
371  // u = u0 + uz * z
372  // v = v0 + vz * z
373  // where u0, uz, v0, vz are the parameters to be obtained from the fit.
374  // These are contained in a vector p which is a 4x1 dim matrix, and
375  // its transpose pT is (u0, v0, uz, vz). Note the ordering!
376  // The covariance matrix for each pair of measurements is 2 x 2 and
377  // the inverse of this is the error matrix E.
378  // The error matrix for the whole set of n measurements is a diagonal
379  // matrix with diagonal elements the individual 2 x 2 error matrices
380  // (because the inverse of a diagonal matrix is a diagonal matrix
381  // with each element the inverse of the original.)
382 
383  // The matrix 'matrix' in method 'CSCSegment::weightMatrix()' is this
384  // block-diagonal overall covariance matrix. It is inverted to the
385  // block-diagonal error matrix right before it is returned.
386 
387  // Use the matrix A defined by
388  // 1 0 z1 0
389  // 0 1 0 z1
390  // 1 0 z2 0
391  // 0 1 0 z2
392  // .. .. .. ..
393  // 1 0 zn 0
394  // 0 1 0 zn
395 
396  // The matrix A is returned by 'CSCSegment::derivativeMatrix()'.
397 
398  // Then the normal equations are encapsulated in the matrix equation
399  //
400  // (AT E A)p = (AT E)m
401  //
402  // where AT is the transpose of A.
403  // We'll call the combined matrix on the LHS, M, and that on the RHS, B:
404  // M p = B
405 
406  // We solve this for the parameter vector, p.
407  // The elements of M and B then involve sums over the hits
408 
409  // The 4 values in p are returned by 'CSCSegment::parameters()'
410  // in the order uz, vz, u0, v0.
411  // The error matrix of the parameters is obtained by
412  // (AT E A)^-1
413  // calculated in 'CSCSegment::parametersErrors()'.
414 
415  // NOTE 1
416  // It does the #hits = 2 case separately.
417  // (I hope they're not on the same layer! They should not be, by construction.)
418 
419  // NOTE 2
420  // We need local position of a RecHit w.r.t. the CHAMBER
421  // and the RecHit itself only knows its local position w.r.t.
422  // the LAYER, so we must explicitly transform global position.
423 
424  CLHEP::HepMatrix M(4,4,0);
425  CLHEP::HepVector B(4,0);
426 
427  ChamberHitContainer::const_iterator ih = proto_segment.begin();
428 
429  for (ih = proto_segment.begin(); ih != proto_segment.end(); ++ih) {
430 
431  const CSCRecHit2D& hit = (**ih);
432  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
433  GlobalPoint gp = layer->toGlobal(hit.localPosition());
434  LocalPoint lp = theChamber->toLocal(gp); // FIX !!
435 
436  // ptc: Local position of hit w.r.t. chamber
437  double u = lp.x();
438  double v = lp.y();
439  double z = lp.z();
440 
441  // ptc: Covariance matrix of local errors MUST BE CHECKED IF COMAPTIBLE
442  CLHEP::HepMatrix IC(2,2);
443  IC(1,1) = hit.localPositionError().xx();
444  IC(1,2) = hit.localPositionError().xy();
445  IC(2,1) = IC(1,2); // since Cov is symmetric
446  IC(2,2) = hit.localPositionError().yy();
447 
448  // ptc: Invert covariance matrix (and trap if it fails!)
449  int ierr = 0;
450  IC.invert(ierr); // inverts in place
451  if (ierr != 0) {
452  LogDebug("CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC << "\n";
453 
454  // @@ NOW WHAT TO DO? Exception? Return? Ignore?
455  }
456 
457  // ptc: Note that IC is no longer Cov but Cov^-1
458  M(1,1) += IC(1,1);
459  M(1,2) += IC(1,2);
460  M(1,3) += IC(1,1) * z;
461  M(1,4) += IC(1,2) * z;
462  B(1) += u * IC(1,1) + v * IC(1,2);
463 
464  M(2,1) += IC(2,1);
465  M(2,2) += IC(2,2);
466  M(2,3) += IC(2,1) * z;
467  M(2,4) += IC(2,2) * z;
468  B(2) += u * IC(2,1) + v * IC(2,2);
469 
470  M(3,1) += IC(1,1) * z;
471  M(3,2) += IC(1,2) * z;
472  M(3,3) += IC(1,1) * z * z;
473  M(3,4) += IC(1,2) * z * z;
474  B(3) += ( u * IC(1,1) + v * IC(1,2) ) * z;
475 
476  M(4,1) += IC(2,1) * z;
477  M(4,2) += IC(2,2) * z;
478  M(4,3) += IC(2,1) * z * z;
479  M(4,4) += IC(2,2) * z * z;
480  B(4) += ( u * IC(2,1) + v * IC(2,2) ) * z;
481  }
482 
483  // Solve the matrix equation using CLHEP's 'solve'
484  //@@ ptc: CAN solve FAIL?? UNCLEAR FROM (LACK OF) CLHEP DOC
485  CLHEP::HepVector p = solve(M, B);
486 
487  // Update member variables uz, vz, theOrigin
488  theOrigin = LocalPoint(p(1), p(2), 0.);
489  uz = p(3);
490  vz = p(4);
491 }
492 
494 
495  // The chi-squared is (m-Ap)T E (m-Ap)
496  // where T denotes transpose.
497  // This collapses to a simple sum over contributions from each
498  // pair of measurements.
499  float u0 = theOrigin.x();
500  float v0 = theOrigin.y();
501  double chsq = 0.;
502 
503  ChamberHitContainer::const_iterator ih;
504  for (ih = proto_segment.begin(); ih != proto_segment.end(); ++ih) {
505 
506  const CSCRecHit2D& hit = (**ih);
507  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
508  GlobalPoint gp = layer->toGlobal(hit.localPosition());
509  LocalPoint lp = theChamber->toLocal(gp); // FIX !!
510 
511  double hu = lp.x();
512  double hv = lp.y();
513  double hz = lp.z();
514 
515  double du = u0 + uz * hz - hu;
516  double dv = v0 + vz * hz - hv;
517 
518  CLHEP::HepMatrix IC(2,2);
519  IC(1,1) = hit.localPositionError().xx();
520  IC(1,2) = hit.localPositionError().xy();
521  IC(2,1) = IC(1,2);
522  IC(2,2) = hit.localPositionError().yy();
523 
524  // Invert covariance matrix
525  int ierr = 0;
526  IC.invert(ierr);
527  if (ierr != 0) {
528  LogDebug("CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC << "\n";
529 
530  // @@ NOW WHAT TO DO? Exception? Return? Ignore?
531  }
532 
533  chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
534  }
535  theChi2 = chsq;
536 }
537 
539 
540  // Always enforce direction of segment to point from IP outwards
541  // (Incorrect for particles not coming from IP, of course.)
542 
543  double dxdz = uz;
544  double dydz = vz;
545  double dz = 1./sqrt(1. + dxdz*dxdz + dydz*dydz);
546  double dx = dz*dxdz;
547  double dy = dz*dydz;
548  LocalVector localDir(dx,dy,dz);
549 
550  // localDir may need sign flip to ensure it points outward from IP
551  // ptc: Examine its direction and origin in global z: to point outward
552  // the localDir should always have same sign as global z...
553 
554  double globalZpos = ( theChamber->toGlobal( theOrigin ) ).z();
555  double globalZdir = ( theChamber->toGlobal( localDir ) ).z();
556  double directionSign = globalZpos * globalZdir;
557 
558  theDirection = (directionSign * localDir).unit();
559 }
560 
561 float CSCSegAlgoTC::phiAtZ(float z) const {
562 
563  // Returns a phi in [ 0, 2*pi )
564  const CSCLayer* l1 = theChamber->layer(1);
565  GlobalPoint gp = l1->toGlobal(theOrigin);
567 
568  float x = gp.x() + (gv.x()/gv.z())*(z - gp.z());
569  float y = gp.y() + (gv.y()/gv.z())*(z - gp.z());
570  float phi = atan2(y, x);
571  if (phi < 0.f)
572  phi += 2. * M_PI;
573 
574  return phi ;
575 }
576 
578 
579  // compare the chi2 of two segments
580  double oldChi2 = theChi2;
581  LocalPoint oldOrigin = theOrigin;
582  LocalVector oldDirection = theDirection;
583  ChamberHitContainer oldSegment = proto_segment;
584 
585  bool ok = replaceHit(h, layer);
586 
587  if (ok) {
588  LogDebug("CSC") << " hit in same layer as a hit on segment; try replacing old one..."
589  << " chi2 new: " << theChi2 << " old: " << oldChi2 << "\n";
590  }
591 
592  if ((theChi2 < oldChi2) && (ok)) {
593  LogDebug("CSC") << " segment with replaced hit is better.\n";
594  }
595  else {
596  proto_segment = oldSegment;
597  theChi2 = oldChi2;
598  theOrigin = oldOrigin;
599  theDirection = oldDirection;
600  }
601 }
602 
604 
605  double oldChi2 = theChi2;
606  LocalPoint oldOrigin = theOrigin;
607  LocalVector oldDirection = theDirection;
608  ChamberHitContainer oldSegment = proto_segment;
609 
610  bool ok = addHit(h, layer);
611 
612  if (ok) {
613  LogDebug("CSC") << " hit in new layer: added to segment, new chi2: "
614  << theChi2 << "\n";
615  }
616 
617  int ndf = 2*proto_segment.size() - 4;
618 
619  if (ok && ((ndf <= 0) || (theChi2/ndf < chi2Max))) {
620  LogDebug("CSC") << " segment with added hit is good.\n" ;
621  }
622  else {
623  proto_segment = oldSegment;
624  theChi2 = oldChi2;
625  theOrigin = oldOrigin;
626  theDirection = oldDirection;
627  }
628 }
629 
631  float h1x = h1->localPosition().x();
632  float h2x = h2->localPosition().x();
633  float deltaX = (h1->localPosition()-h2->localPosition()).x();
634  LogDebug("CSC") << " Hits at local x= " << h1x << ", "
635  << h2x << " have separation= " << deltaX;
636  return (fabs(deltaX) < (dRPhiMax))? true:false; // +v
637 }
638 
640 
641  const CSCLayer* l1 = theChamber->layer(h1->cscDetId().layer());
642  GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
643  const CSCLayer* l2 = theChamber->layer(h2->cscDetId().layer());
644  GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
645 
646  float h1p = gp1.phi();
647  float h2p = gp2.phi();
648  float dphi12 = h1p - h2p;
649 
650  // Into range [-pi, pi) (phi() returns values in this range)
651  if (dphi12 < -M_PI)
652  dphi12 += 2.*M_PI;
653  if (dphi12 > M_PI)
654  dphi12 -= 2.*M_PI;
655  LogDebug("CSC") << " Hits at global phi= " << h1p << ", "
656  << h2p << " have separation= " << dphi12;
657  return (fabs(dphi12) < dPhiMax)? true:false; // +v
658 }
659 
661 
662  // Is hit near segment?
663  // Requires deltaphi and rxy*deltaphi within ranges specified
664  // in orcarc, or by default, where rxy=sqrt(x**2+y**2) of hit itself.
665  // Note that to make intuitive cuts on delta(phi) one must work in
666  // phi range (-pi, +pi] not [0, 2pi
667 
668  const CSCLayer* l1 = theChamber->layer(h->cscDetId().layer());
669  GlobalPoint hp = l1->toGlobal(h->localPosition());
670 
671  float hphi = hp.phi(); // in (-pi, +pi]
672  if (hphi < 0.)
673  hphi += 2.*M_PI; // into range [0, 2pi)
674  float sphi = phiAtZ(hp.z()); // in [0, 2*pi)
675  float phidif = sphi-hphi;
676  if (phidif < 0.)
677  phidif += 2.*M_PI; // into range [0, 2pi)
678  if (phidif > M_PI)
679  phidif -= 2.*M_PI; // into range (-pi, pi]
680 
681  float dRPhi = fabs(phidif)*hp.perp();
682  LogDebug("CSC") << " is hit at phi_h= " << hphi << " near segment phi_seg= " << sphi
683  << "? is " << dRPhi << "<" << dRPhiFineMax << " ? "
684  << " and is |" << phidif << "|<" << dPhiFineMax << " ?";
685 
686  return ((dRPhi < dRPhiFineMax) &&
687  (fabs(phidif) < dPhiFineMax))? true:false; // +v
688 }
689 
690 bool CSCSegAlgoTC::hasHitOnLayer(int layer) const {
691 
692  // Is there is already a hit on this layer?
694 
695  for(it = proto_segment.begin(); it != proto_segment.end(); it++)
696  if ((*it)->cscDetId().layer() == layer)
697  return true;
698 
699  return false;
700 }
701 
703 
704  // Dump positions of RecHit's in each CSCChamber
706 
707  for(it=rechits.begin(); it!=rechits.end(); it++) {
708 
709  const CSCLayer* l1 = theChamber->layer((*it)->cscDetId().layer());
710  GlobalPoint gp1 = l1->toGlobal((*it)->localPosition());
711 
712  LogDebug("CSC") << "Global pos.: " << gp1 << ", phi: " << gp1.phi() << ". Local position: "
713  << (*it)->localPosition() << ", phi: "
714  << (*it)->localPosition().phi() << ". Layer: "
715  << (*it)->cscDetId().layer() << "\n";
716  }
717 }
718 
719 
720 bool CSCSegAlgoTC::isSegmentGood(std::vector<ChamberHitContainer>::iterator seg, std::vector<double>::iterator chi2,
721  const ChamberHitContainer& rechitsInChamber, BoolContainer& used) const {
722 
723  // Apply any selection cuts to segment
724 
725  // 1) Require a minimum no. of hits
726  // (@@ THESE VALUES SHOULD BECOME PARAMETERS?)
727 
728  // 2) Ensure no hits on segment are already assigned to another segment
729  // (typically of higher quality)
730 
731  unsigned int iadd = (rechitsInChamber.size() > 20 )? 1 : 0;
732 
733  if (seg->size() < 3 + iadd)
734  return false;
735 
736  // Additional part of alternative segment selection scheme: reject
737  // segments with a chi2 probability of less than chi2ndfProbMin. Relies on list
738  // being sorted with "SegmentSorting == 2", that is first by nrechits and then
739  // by chi2prob in subgroups of same nr of rechits.
740 
741  if( SegmentSorting == 2 ){
742  if( (*chi2) != 0 && ((2*seg->size())-4) >0 ) {
743  if ( ChiSquaredProbability((*chi2),(double)(2*seg->size()-4)) < chi2ndfProbMin ) {
744  return false;
745  }
746  }
747  if((*chi2) == 0 ) return false;
748  }
749 
750 
751  for(unsigned int ish = 0; ish < seg->size(); ++ish) {
752 
753  ChamberHitContainerCIt ib = rechitsInChamber.begin();
754 
755  for(ChamberHitContainerCIt ic = ib; ic != rechitsInChamber.end(); ++ic) {
756 
757  if(((*seg)[ish] == (*ic)) && used[ic-ib])
758  return false;
759  }
760  }
761 
762  return true;
763 }
764 
765 void CSCSegAlgoTC::flagHitsAsUsed(std::vector<ChamberHitContainer>::iterator seg, const ChamberHitContainer& rechitsInChamber,
766  BoolContainer& used) const {
767 
768  // Flag hits on segment as used
769 
770  ChamberHitContainerCIt ib = rechitsInChamber.begin();
771 
772  for(unsigned int ish = 0; ish < seg->size(); ish++) {
773 
774  for(ChamberHitContainerCIt iu = ib; iu != rechitsInChamber.end(); ++iu)
775  if((*seg)[ish] == (*iu))
776  used[iu-ib] = true;
777  }
778 }
779 
781 
782  // Sort the segment store according to segment 'quality' (chi2/#hits ?) and
783  // remove any segments which contain hits assigned to higher-quality segments.
784 
785  if (candidates.empty())
786  return;
787 
788  // Initialize flags that a given hit has been allocated to a segment
789  BoolContainer used(rechitsInChamber.size(), false);
790 
791  // Sort by chi2/#hits
792  segmentSort();
793 
794  // Select best quality segments, requiring hits are assigned to just one segment
795  // Because I want to erase the bad segments, the iterator must be incremented
796  // inside the loop, and only when the erase is not called
797 
798  std::vector<ChamberHitContainer>::iterator is;
799  std::vector<double>::iterator ichi = chi2s.begin();
800  std::vector<AlgebraicSymMatrix>::iterator iErr = errors.begin();
801  std::vector<LocalPoint>::iterator iOrig = origins.begin();
802  std::vector<LocalVector>::iterator iDir = directions.begin();
803 
804  for (is = candidates.begin(); is != candidates.end(); ) {
805 
806  bool goodSegment = isSegmentGood(is, ichi, rechitsInChamber, used);
807 
808  if (goodSegment) {
809  LogDebug("CSC") << "Accepting segment: ";
810 
811  flagHitsAsUsed(is, rechitsInChamber, used);
812  ++is;
813  ++ichi;
814  ++iErr;
815  ++iOrig;
816  ++iDir;
817  }
818  else {
819  LogDebug("CSC") << "Rejecting segment: ";
820  is = candidates.erase(is);
821  ichi = chi2s.erase(ichi);
822  iErr = errors.erase(iErr);
823  iOrig = origins.erase(iOrig);
824  iDir = directions.erase(iDir);
825  }
826  }
827 }
828 
830 
831  // The segment collection is sorted according chi2/#hits
832 
833  for(unsigned int i=0; i<candidates.size()-1; i++) {
834  for(unsigned int j=i+1; j<candidates.size(); j++) {
835 
838  if (i == j)
839  continue;
840 
841  int n1 = candidates[i].size();
842  int n2 = candidates[j].size();
843 
844  if( SegmentSorting == 2 ){ // Sort criterion: first sort by Nr of rechits, then in groups of rechits by chi2prob:
845  if ( n2 > n1 ) { // sort by nr of rechits
847  candidates[j] = candidates[i];
848  candidates[i] = temp;
849 
850  double temp1 = chi2s[j];
851  chi2s[j] = chi2s[i];
852  chi2s[i] = temp1;
853 
854  AlgebraicSymMatrix temp2 = errors[j];
855  errors[j] = errors[i];
856  errors[i] = temp2;
857 
858  LocalPoint temp3 = origins[j];
859  origins[j] = origins[i];
860  origins[i] = temp3;
861 
862  LocalVector temp4 = directions[j];
863  directions[j] = directions[i];
864  directions[i] = temp4;
865  }
866  // sort by chi2 probability in subgroups with equal nr of rechits
867  if(chi2s[i] != 0. && 2*n2-4 > 0 ) {
868  if( n2 == n1 && (ChiSquaredProbability( chi2s[i],(double)(2*n1-4)) < ChiSquaredProbability(chi2s[j],(double)(2*n2-4))) ){
870  candidates[j] = candidates[i];
871  candidates[i] = temp;
872 
873  double temp1 = chi2s[j];
874  chi2s[j] = chi2s[i];
875  chi2s[i] = temp1;
876 
877  AlgebraicSymMatrix temp2 = errors[j];
878  errors[j] = errors[i];
879  errors[i] = temp2;
880 
881  LocalPoint temp3 = origins[j];
882  origins[j] = origins[i];
883  origins[i] = temp3;
884 
885  LocalVector temp4 = directions[j];
886  directions[j] = directions[i];
887  directions[i] = temp4;
888  }
889  }
890  }
891  else if( SegmentSorting == 1 ){
892  if ((chi2s[i]/n1) > (chi2s[j]/n2)) {
893 
895  candidates[j] = candidates[i];
896  candidates[i] = temp;
897 
898  double temp1 = chi2s[j];
899  chi2s[j] = chi2s[i];
900  chi2s[i] = temp1;
901 
902  AlgebraicSymMatrix temp2 = errors[j];
903  errors[j] = errors[i];
904  errors[i] = temp2;
905 
906  LocalPoint temp3 = origins[j];
907  origins[j] = origins[i];
908  origins[i] = temp3;
909 
910  LocalVector temp4 = directions[j];
911  directions[j] = directions[i];
912  directions[i] = temp4;
913  }
914  }
915  else{
916  LogDebug("CSC") << "No valid segment sorting specified - BAD !!!\n";
917  }
918  }
919  }
920 }
921 
923 
926 
927  // (AT W A)^-1
928  // from http://www.phys.ufl.edu/~avery/fitting.html, part I
929  int ierr;
930  AlgebraicSymMatrix result = weights.similarityT(A);
931  result.invert(ierr);
932 
933  // blithely assuming the inverting never fails...
934  return result;
935 }
936 
937 CLHEP::HepMatrix CSCSegAlgoTC::derivativeMatrix() const {
938 
939  ChamberHitContainer::const_iterator it;
940  int nhits = proto_segment.size();
941  CLHEP::HepMatrix matrix(2*nhits, 4);
942  int row = 0;
943 
944  for(it = proto_segment.begin(); it != proto_segment.end(); ++it) {
945 
946  const CSCRecHit2D& hit = (**it);
947  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
948  GlobalPoint gp = layer->toGlobal(hit.localPosition());
949  LocalPoint lp = theChamber->toLocal(gp); // FIX
950  float z = lp.z();
951  ++row;
952  matrix(row, 1) = 1.;
953  matrix(row, 3) = z;
954  ++row;
955  matrix(row, 2) = 1.;
956  matrix(row, 4) = z;
957  }
958  return matrix;
959 }
960 
961 
963 
964  std::vector<const CSCRecHit2D*>::const_iterator it;
965  int nhits = proto_segment.size();
966  AlgebraicSymMatrix matrix(2*nhits, 0);
967  int row = 0;
968 
969  for (it = proto_segment.begin(); it != proto_segment.end(); ++it) {
970 
971  const CSCRecHit2D& hit = (**it);
972  ++row;
973  matrix(row, row) = hit.localPositionError().xx();
974  matrix(row, row+1) = hit.localPositionError().xy();
975  ++row;
976  matrix(row, row-1) = hit.localPositionError().xy();
977  matrix(row, row) = hit.localPositionError().yy();
978  }
979  int ierr;
980  matrix.invert(ierr);
981  return matrix;
982 }
983 
985 
986  // The CSCSegment needs the error matrix re-arranged
987 
988  AlgebraicSymMatrix hold( a );
989 
990  // errors on slopes into upper left
991  a(1,1) = hold(3,3);
992  a(1,2) = hold(3,4);
993  a(2,1) = hold(4,3);
994  a(2,2) = hold(4,4);
995 
996  // errors on positions into lower right
997  a(3,3) = hold(1,1);
998  a(3,4) = hold(1,2);
999  a(4,3) = hold(2,1);
1000  a(4,4) = hold(2,2);
1001 
1002  // off-diagonal elements remain unchanged
1003 
1004 }
1005 
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void fillLocalDirection()
int i
Definition: DBlmapReader.cc:9
float xx() const
Definition: LocalError.h:24
CSCDetId cscDetId() const
Definition: CSCRecHit2D.h:52
int ib
Definition: cuy.py:660
T perp() const
Definition: PV3DBase.h:72
const std::string myName
Definition: CSCSegAlgoTC.h:210
double_binary B
Definition: DDStreamer.cc:234
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
std::vector< AlgebraicSymMatrix > errors
Definition: CSCSegAlgoTC.h:159
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:67
susybsm::HSCParticleRefVector hv
Definition: classes.h:28
void compareProtoSegment(const CSCRecHit2D *h, int layer)
AlgebraicSymMatrix weightMatrix() const
std::deque< bool > BoolContainer
Definition: CSCSegAlgoTC.h:57
int layer() const
Definition: CSCDetId.h:74
float phiAtZ(float z) const
AlgebraicSymMatrix calculateError() const
float float float z
std::vector< CSCSegment > buildSegments(const ChamberHitContainer &rechits)
Definition: CSCSegAlgoTC.cc:58
tuple s2
Definition: indexGen.py:106
bool areHitsCloseInLocalX(const CSCRecHit2D *h1, const CSCRecHit2D *h2) const
bool isHitNearSegment(const CSCRecHit2D *h) const
std::string chamberTypeName() const
void updateParameters()
float xy() const
Definition: LocalError.h:25
void dumpHits(const ChamberHitContainer &rechits) const
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:46
const CSCChamberSpecs * specs() const
Definition: CSCChamber.h:42
CLHEP::HepMatrix AlgebraicMatrix
std::vector< double > chi2s
Definition: CSCSegAlgoTC.h:160
LocalError localPositionError() const
Definition: CSCRecHit2D.h:51
float yy() const
Definition: LocalError.h:26
string unit
Definition: csvLumiCalc.py:46
susybsm::HSCParticleRefProd hp
Definition: classes.h:27
T sqrt(T t)
Definition: SSEVec.h:48
T phi() const
Definition: Phi.h:41
LocalVector theDirection
Definition: CSCSegAlgoTC.h:165
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
void fillChiSquared()
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:39
void flagHitsAsUsed(std::vector< ChamberHitContainer >::iterator is, const ChamberHitContainer &rechitsInChamber, BoolContainer &used) const
float dRPhiFineMax
Definition: CSCSegAlgoTC.h:180
double f[11][100]
float ChiSquaredProbability(double chiSquared, double nrDOF)
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
bool isSegmentGood(std::vector< ChamberHitContainer >::iterator is, std::vector< double >::iterator ichi, const ChamberHitContainer &rechitsInChamber, BoolContainer &used) const
bool areHitsCloseInGlobalPhi(const CSCRecHit2D *h1, const CSCRecHit2D *h2) const
#define M_PI
bool addHit(const CSCRecHit2D *aHit, int layer)
Utility functions.
float chi2ndfProbMin
Definition: CSCSegAlgoTC.h:175
void segmentSort()
double theChi2
Definition: CSCSegAlgoTC.h:163
const CSCChamber * theChamber
Member variables.
Definition: CSCSegAlgoTC.h:155
std::vector< CSCSegment > run(const CSCChamber *aChamber, const ChamberHitContainer &rechits)
Definition: CSCSegAlgoTC.cc:53
ChamberHitContainer proto_segment
Definition: CSCSegAlgoTC.h:162
std::vector< const CSCRecHit2D * > ChamberHitContainer
Definition: CSCSegAlgoTC.h:49
std::vector< LocalPoint > origins
Definition: CSCSegAlgoTC.h:157
void tryAddingHitsToSegment(const ChamberHitContainer &rechits, const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2)
void fitSlopes()
bool hasHitOnLayer(int layer) const
CLHEP::HepMatrix derivativeMatrix() const
LocalPoint theOrigin
Definition: CSCSegAlgoTC.h:164
double a
Definition: hdecay.h:121
CLHEP::HepSymMatrix AlgebraicSymMatrix
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
std::vector< int > LayerIndex
Typedefs.
Definition: CSCSegAlgoTC.h:48
std::vector< ChamberHitContainer > candidates
Definition: CSCSegAlgoTC.h:156
Definition: DDAxes.h:10
LocalPoint localPosition() const
Definition: CSCRecHit2D.h:50
T x() const
Definition: PV3DBase.h:62
ChamberHitContainer::const_iterator ChamberHitContainerCIt
Definition: CSCSegAlgoTC.h:50
void increaseProtoSegment(const CSCRecHit2D *h, int layer)
void pruneTheSegments(const ChamberHitContainer &rechitsInChamber)
void flipErrors(AlgebraicSymMatrix &) const
bool replaceHit(const CSCRecHit2D *h, int layer)
CSCSegAlgoTC(const edm::ParameterSet &ps)
Constructor.
Definition: CSCSegAlgoTC.cc:27
float dPhiFineMax
Definition: CSCSegAlgoTC.h:185
Definition: DDAxes.h:10