CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCSegAlgoSK.cc
Go to the documentation of this file.
1 
9 #include "CSCSegAlgoSK.h"
11 // For clhep Matrix::solve
14 
17 
18 #include <algorithm>
19 #include <cmath>
20 #include <iostream>
21 #include <string>
22 
24  myName("CSCSegAlgoSK") {
25 
26  debugInfo = ps.getUntrackedParameter<bool>("verboseInfo");
27 
28  dRPhiMax = ps.getParameter<double>("dRPhiMax");
29  dPhiMax = ps.getParameter<double>("dPhiMax");
30  dRPhiFineMax = ps.getParameter<double>("dRPhiFineMax");
31  dPhiFineMax = ps.getParameter<double>("dPhiFineMax");
32  chi2Max = ps.getParameter<double>("chi2Max");
33  wideSeg = ps.getParameter<double>("wideSeg");
34  minLayersApart = ps.getParameter<int>("minLayersApart");
35 
36  LogDebug("CSC") << myName << " has algorithm cuts set to: \n"
37  << "--------------------------------------------------------------------\n"
38  << "dRPhiMax = " << dRPhiMax << '\n'
39  << "dPhiMax = " << dPhiMax << '\n'
40  << "dRPhiFineMax = " << dRPhiFineMax << '\n'
41  << "dPhiFineMax = " << dPhiFineMax << '\n'
42  << "chi2Max = " << chi2Max << '\n'
43  << "wideSeg = " << wideSeg << '\n'
44  << "minLayersApart = " << minLayersApart << std::endl;
45 }
46 
47 std::vector<CSCSegment> CSCSegAlgoSK::run(const CSCChamber* aChamber, ChamberHitContainer rechits) {
48  theChamber = aChamber;
49  return buildSegments(rechits);
50 }
51 
53 
54  LogDebug("CSC") << "*********************************************";
55  LogDebug("CSC") << "Start segment building in the new chamber: " << theChamber->specs()->chamberTypeName();
56  LogDebug("CSC") << "*********************************************";
57 
58 
59  LayerIndex layerIndex(rechits.size());
60 
61  for(unsigned int i = 0; i < rechits.size(); i++) {
62 
63  layerIndex[i] = rechits[i]->cscDetId().layer();
64  }
65 
66  double z1 = theChamber->layer(1)->position().z();
67  double z6 = theChamber->layer(6)->position().z();
68 
69  if ( z1 > 0. ) {
70  if ( z1 > z6 ) {
71  reverse(layerIndex.begin(), layerIndex.end());
72  reverse(rechits.begin(), rechits.end());
73  }
74  }
75  else if ( z1 < 0. ) {
76  if ( z1 < z6 ) {
77  reverse(layerIndex.begin(), layerIndex.end());
78  reverse(rechits.begin(), rechits.end());
79  }
80  }
81 
82  if (debugInfo) {
83  // dump after sorting
84  dumpHits(rechits);
85  }
86 
87  if (rechits.size() < 2) {
88  LogDebug("CSC") << myName << ": " << rechits.size() <<
89  " hit(s) in chamber is not enough to build a segment.\n";
90  return std::vector<CSCSegment>();
91  }
92 
93  // We have at least 2 hits. We intend to try all possible pairs of hits to start
94  // segment building. 'All possible' means each hit lies on different layers in the chamber.
95  // BUT... once a hit has been assigned to a segment, we don't consider
96  // it again.
97 
98  // Choose first hit (as close to IP as possible) h1 and
99  // second hit (as far from IP as possible) h2
100  // To do this we iterate over hits in the chamber by layer - pick two layers.
101  // @@ Require the two layers are at least 3 layers apart. May need tuning?
102  // Then we iterate over hits within each of these layers and pick h1 and h2 from these.
103  // If they are 'close enough' we build an empty segment.
104  // Then try adding hits to this segment.
105 
106  // Initialize flags that a given hit has been allocated to a segment
107  BoolContainer used(rechits.size(), false);
108 
109  // Define buffer for segments we build
110  std::vector<CSCSegment> segments;
111 
112  ChamberHitContainerCIt ib = rechits.begin();
113  ChamberHitContainerCIt ie = rechits.end();
114 
115  // Possibly allow 2 passes, second widening scale factor for cuts
116  windowScale = 1.; // scale factor for cuts
117 
118  int npass = (wideSeg > 1.)? 2 : 1;
119 
120  for (int ipass = 0; ipass < npass; ++ipass) {
121  for (ChamberHitContainerCIt i1 = ib; i1 != ie; ++i1) {
122  bool segok = false;
123  if(used[i1-ib])
124  continue;
125 
126  int layer1 = layerIndex[i1-ib]; //(*i1)->cscDetId().layer();
127  const CSCRecHit2D* h1 = *i1;
128 
129  for (ChamberHitContainerCIt i2 = ie-1; i2 != i1; --i2) {
130  if(used[i2-ib])
131  continue;
132 
133  int layer2 = layerIndex[i2-ib]; //(*i2)->cscDetId().layer();
134 
135  if (abs(layer2 - layer1) < minLayersApart)
136  break;
137  const CSCRecHit2D* h2 = *i2;
138 
139  if (areHitsCloseInLocalX(h1, h2) && areHitsCloseInGlobalPhi(h1, h2)) {
140 
141  proto_segment.clear();
142 
143  const CSCLayer* l1 = theChamber->layer(layer1);
144  GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
145  const CSCLayer* l2 = theChamber->layer(layer2);
146  GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
147  LogDebug("CSC") << "start new segment from hits " << "h1: "
148  << gp1 << " - h2: " << gp2 << "\n";
149 
150  if (!addHit(h1, layer1)) {
151  LogDebug("CSC") << " failed to add hit h1\n";
152  continue;
153  }
154 
155  if (!addHit(h2, layer2)) {
156  LogDebug("CSC") << " failed to add hit h2\n";
157  continue;
158  }
159 
160  tryAddingHitsToSegment(rechits, used, layerIndex, i1, i2);
161 
162  // Check no. of hits on segment, and if enough flag them as used
163  // and store the segment
164  segok = isSegmentGood(rechits);
165  if (segok) {
166  flagHitsAsUsed(rechits, used);
167  // Copy the proto_segment and its properties inside a CSCSegment.
168  // Then fill the segment vector..
169 
170  if (proto_segment.empty()) {
171  LogDebug("CSC") << "No segment has been found !!!\n";
172  }
173  else {
174 
175  // calculate error matrix...
177 
178  // but reorder components to match what's required by TrackingRecHit interface
179  // i.e. slopes first, then positions
180 
181  flipErrors( errors );
182 
184 
185  LogDebug("CSC") << "Found a segment !!!\n";
186  segments.push_back(temp);
187  }
188  }
189  } // h1 & h2 close
190 
191  if (segok)
192  break;
193  } // i2
194  } // i1
195 
196  if (segments.size() > 1)
197  break; // only change window if no segments found
198 
199  // Increase cut windows by factor of wideSeg
201 
202  } // ipass
203 
204  // Give the segments to the CSCChamber
205  return segments;
206 }
207 
209  BoolContainer used, LayerIndex layerIndex,
210  const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2) {
211 
212  // Iterate over the layers with hits in the chamber
213  // Skip the layers containing the segment endpoints
214  // Test each hit on the other layers to see if it is near the segment
215  // If it is, see whether there is already a hit on the segment from the same layer
216  // - if so, and there are more than 2 hits on the segment, copy the segment,
217  // replace the old hit with the new hit. If the new segment chi2 is better
218  // then replace the original segment with the new one (by swap)
219  // - if not, copy the segment, add the hit. If the new chi2/dof is still satisfactory
220  // then replace the original segment with the new one (by swap)
221 
222  ChamberHitContainerCIt ib = rechits.begin();
223  ChamberHitContainerCIt ie = rechits.end();
224 
225  for (ChamberHitContainerCIt i = ib; i != ie; ++i) {
226  if (i == i1 || i == i2 || used[i-ib])
227  continue;
228 
229  int layer = layerIndex[i-ib];
230  const CSCRecHit2D* h = *i;
231  if (isHitNearSegment(h)) {
232 
233  GlobalPoint gp1 = theChamber->layer(layer)->toGlobal(h->localPosition());
234  LogDebug("CSC") << " hit at global " << gp1 << " is near segment\n.";
235 
236  // Don't consider alternate hits on layers holding the two starting points
237  if (hasHitOnLayer(layer)) {
238  if (proto_segment.size() <= 2) {
239  LogDebug("CSC") << " " << proto_segment.size()
240  << " hits on segment...skip hit on same layer.\n";
241  continue;
242  }
243  compareProtoSegment(h, layer);
244  }
245  else
246  increaseProtoSegment(h, layer);
247  } // h & seg close
248  } // i
249 }
250 
252  float h1x = h1->localPosition().x();
253  float h2x = h2->localPosition().x();
254  float deltaX = (h1->localPosition()-h2->localPosition()).x();
255  LogDebug("CSC") << " Hits at local x= " << h1x << ", "
256  << h2x << " have separation= " << deltaX;
257  return (fabs(deltaX) < (dRPhiMax * windowScale))? true:false; // +v
258 }
259 
261 
262  const CSCLayer* l1 = theChamber->layer(h1->cscDetId().layer());
263  GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
264  const CSCLayer* l2 = theChamber->layer(h2->cscDetId().layer());
265  GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
266 
267  float h1p = gp1.phi();
268  float h2p = gp2.phi();
269  float dphi12 = h1p - h2p;
270 
271  // Into range [-pi, pi) (phi() returns values in this range)
272  if (dphi12 < -M_PI)
273  dphi12 += 2.*M_PI;
274  if (dphi12 > M_PI)
275  dphi12 -= 2.*M_PI;
276  LogDebug("CSC") << " Hits at global phi= " << h1p << ", "
277  << h2p << " have separation= " << dphi12;
278  return (fabs(dphi12) < (dPhiMax * windowScale))? true:false; // +v
279 }
280 
282 
283  // Is hit near segment?
284  // Requires deltaphi and rxy*deltaphi within ranges specified
285  // in parameter set, where rxy=sqrt(x**2+y**2) of hit itself.
286  // Note that to make intuitive cuts on delta(phi) one must work in
287  // phi range (-pi, +pi] not [0, 2pi)
288 
289  const CSCLayer* l1 = theChamber->layer(h->cscDetId().layer());
290  GlobalPoint hp = l1->toGlobal(h->localPosition());
291 
292  float hphi = hp.phi(); // in (-pi, +pi]
293  if (hphi < 0.)
294  hphi += 2.*M_PI; // into range [0, 2pi)
295  float sphi = phiAtZ(hp.z()); // in [0, 2*pi)
296  float phidif = sphi-hphi;
297  if (phidif < 0.)
298  phidif += 2.*M_PI; // into range [0, 2pi)
299  if (phidif > M_PI)
300  phidif -= 2.*M_PI; // into range (-pi, pi]
301 
302  float dRPhi = fabs(phidif)*hp.perp();
303  LogDebug("CSC") << " is hit at phi_h= " << hphi << " near segment phi_seg= " << sphi
304  << "? is " << dRPhi << "<" << dRPhiFineMax*windowScale << " ? "
305  << " and is |" << phidif << "|<" << dPhiFineMax*windowScale << " ?";
306 
307  return ((dRPhi < dRPhiFineMax*windowScale) &&
308  (fabs(phidif) < dPhiFineMax*windowScale))? true:false; // +v
309 }
310 
311 float CSCSegAlgoSK::phiAtZ(float z) const {
312 
313  // Returns a phi in [ 0, 2*pi )
314  const CSCLayer* l1 = theChamber->layer((*(proto_segment.begin()))->cscDetId().layer());
315  GlobalPoint gp = l1->toGlobal(theOrigin);
317 
318  float x = gp.x() + (gv.x()/gv.z())*(z - gp.z());
319  float y = gp.y() + (gv.y()/gv.z())*(z - gp.z());
320  float phi = atan2(y, x);
321  if (phi < 0.f )
322  phi += 2. * M_PI;
323 
324  return phi ;
325 }
326 
328 
329  // Dump positions of RecHit's in each CSCChamber
331  edm::LogInfo("CSC") << "CSCChamber rechit dump.\n";
332  for(it=rechits.begin(); it!=rechits.end(); it++) {
333 
334  const CSCLayer* l1 = theChamber->layer((*it)->cscDetId().layer());
335  GlobalPoint gp1 = l1->toGlobal((*it)->localPosition());
336 
337  edm::LogInfo("CSC") << "Global pos.: " << gp1 << ", phi: " << gp1.phi() << ". Local position: "
338  << (*it)->localPosition() << ", phi: "
339  << (*it)->localPosition().phi() << ". Layer: "
340  << (*it)->cscDetId().layer() << "\n";
341  }
342 }
343 
344 bool CSCSegAlgoSK::isSegmentGood(const ChamberHitContainer& rechitsInChamber) const {
345 
346  // If the chamber has 20 hits or fewer, require at least 3 hits on segment
347  // If the chamber has >20 hits require at least 4 hits
348  //@@ THESE VALUES SHOULD BECOME PARAMETERS?
349  bool ok = false;
350 
351  unsigned int iadd = ( rechitsInChamber.size() > 20 )? 1 : 0;
352 
353  if (windowScale > 1.)
354  iadd = 1;
355 
356  if (proto_segment.size() >= 3+iadd)
357  ok = true;
358 
359  return ok;
360 }
361 
362 void CSCSegAlgoSK::flagHitsAsUsed(const ChamberHitContainer& rechitsInChamber,
363  BoolContainer& used ) const {
364 
365  // Flag hits on segment as used
366  ChamberHitContainerCIt ib = rechitsInChamber.begin();
367  ChamberHitContainerCIt hi, iu;
368 
369  for(hi = proto_segment.begin(); hi != proto_segment.end(); ++hi) {
370  for(iu = ib; iu != rechitsInChamber.end(); ++iu) {
371  if(*hi == *iu)
372  used[iu-ib] = true;
373  }
374  }
375 }
376 
377 bool CSCSegAlgoSK::addHit(const CSCRecHit2D* aHit, int layer) {
378 
379  // Return true if hit was added successfully
380  // (and then parameters are updated).
381  // Return false if there is already a hit on the same layer, or insert failed.
382 
383  ChamberHitContainer::const_iterator it;
384 
385  for(it = proto_segment.begin(); it != proto_segment.end(); it++)
386  if (((*it)->cscDetId().layer() == layer) && (aHit != (*it)))
387  return false;
388 
389  proto_segment.push_back(aHit);
391 
392  return true;
393 }
394 
396 
397  // Note that we need local position of a RecHit w.r.t. the CHAMBER
398  // and the RecHit itself only knows its local position w.r.t.
399  // the LAYER, so need to explicitly transform to global position.
400 
401  // no. of hits in the RecHitsOnSegment
402  // By construction this is the no. of layers with hits
403  // since we allow just one hit per layer in a segment.
404 
405  int nh = proto_segment.size();
406 
407  // First hit added to a segment must always fail here
408  if (nh < 2)
409  return;
410 
411  if (nh == 2) {
412 
413  // Once we have two hits we can calculate a straight line
414  // (or rather, a straight line for each projection in xz and yz.)
415  ChamberHitContainer::const_iterator ih = proto_segment.begin();
416  int il1 = (*ih)->cscDetId().layer();
417  const CSCRecHit2D& h1 = (**ih);
418  ++ih;
419  int il2 = (*ih)->cscDetId().layer();
420  const CSCRecHit2D& h2 = (**ih);
421 
422  //@@ Skip if on same layer, but should be impossible
423  if (il1 == il2)
424  return;
425 
426  const CSCLayer* layer1 = theChamber->layer(il1);
427  const CSCLayer* layer2 = theChamber->layer(il2);
428 
429  GlobalPoint h1glopos = layer1->toGlobal(h1.localPosition());
430  GlobalPoint h2glopos = layer2->toGlobal(h2.localPosition());
431 
432  // localPosition is position of hit wrt layer (so local z = 0)
433  theOrigin = h1.localPosition();
434 
435  // We want hit wrt chamber (and local z will be != 0)
436  LocalPoint h1pos = theChamber->toLocal(h1glopos);
437  LocalPoint h2pos = theChamber->toLocal(h2glopos);
438 
439  float dz = h2pos.z()-h1pos.z();
440  uz = (h2pos.x()-h1pos.x())/dz ;
441  vz = (h2pos.y()-h1pos.y())/dz ;
442 
443  theChi2 = 0.;
444  }
445  else if (nh > 2) {
446 
447  // When we have more than two hits then we can fit projections to straight lines
448  fitSlopes();
449  fillChiSquared();
450  } // end of 'if' testing no. of hits
451 
453 }
454 
456 
457  // Update parameters of fit
458  // ptc 13-Aug-02: This does a linear least-squares fit
459  // to the hits associated with the segment, in the z projection.
460 
461  // In principle perhaps one could fit the strip and wire
462  // measurements (u, v respectively), to
463  // u = u0 + uz * z
464  // v = v0 + vz * z
465  // where u0, uz, v0, vz are the parameters resulting from the fit.
466  // But what is actually done is fit to the local x, y coordinates
467  // of the RecHits. However the strip measurement controls the precision
468  // of x, and the wire measurement controls that of y.
469  // Precision in local coordinate:
470  // u (strip, sigma~200um), v (wire, sigma~1cm)
471 
472  // I have verified that this code agrees with the formulation given
473  // on p246-247 of 'Data analysis techniques for high-energy physics
474  // experiments' by Bock, Grote, Notz & Regler, and that on p111-113
475  // of 'Statistics' by Barlow.
476 
477  // Formulate the matrix equation representing the least-squares fit
478  // We have a vector of measurements m, which is a 2n x 1 dim matrix
479  // The transpose mT is (u1, v1, u2, v2, ..., un, vn)
480  // where ui is the strip-associated measurement and vi is the
481  // wire-associated measurement for a given RecHit i.
482  // The fit is to
483  // u = u0 + uz * z
484  // v = v0 + vz * z
485  // where u0, uz, v0, vz are the parameters to be obtained from the fit.
486  // These are contained in a vector p which is a 4x1 dim matrix, and
487  // its transpose pT is (u0, v0, uz, vz). Note the ordering!
488  // The covariance matrix for each pair of measurements is 2 x 2 and
489  // the inverse of this is the error matrix E.
490  // The error matrix for the whole set of n measurements is a diagonal
491  // matrix with diagonal elements the individual 2 x 2 error matrices
492  // (because the inverse of a diagonal matrix is a diagonal matrix
493  // with each element the inverse of the original.)
494 
495  // The matrix 'matrix' in method 'CSCSegment::weightMatrix()' is this
496  // block-diagonal overall covariance matrix. It is inverted to the
497  // block-diagonal error matrix right before it is returned.
498 
499  // Use the matrix A defined by
500  // 1 0 z1 0
501  // 0 1 0 z1
502  // 1 0 z2 0
503  // 0 1 0 z2
504  // .. .. .. ..
505  // 1 0 zn 0
506  // 0 1 0 zn
507 
508  // The matrix A is returned by 'CSCSegment::derivativeMatrix()'.
509 
510  // Then the normal equations are encapsulated in the matrix equation
511  //
512  // (AT E A)p = (AT E)m
513  //
514  // where AT is the transpose of A.
515  // We'll call the combined matrix on the LHS, M, and that on the RHS, B:
516  // M p = B
517 
518  // We solve this for the parameter vector, p.
519  // The elements of M and B then involve sums over the hits
520 
521  // The error matrix of the parameters is obtained by
522  // (AT E A)^-1 calculated in 'calculateError()'.
523 
524  // The 4 values in p can be accessed from 'CSCSegment::parameters()'
525  // in the order uz, vz, u0, v0.
526 
527  // NOTE 1
528  // Do the #hits = 2 case separately.
529  // (I hope they're not on the same layer! They should not be, by construction.)
530 
531  // NOTE 2
532  // We need local position of a RecHit w.r.t. the CHAMBER
533  // and the RecHit itself only knows its local position w.r.t.
534  // the LAYER, so we must explicitly transform global position.
535 
536  CLHEP::HepMatrix M(4,4,0);
537  CLHEP::HepVector B(4,0);
538 
539  ChamberHitContainer::const_iterator ih = proto_segment.begin();
540 
541  for (ih = proto_segment.begin(); ih != proto_segment.end(); ++ih) {
542 
543  const CSCRecHit2D& hit = (**ih);
544  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
545  GlobalPoint gp = layer->toGlobal(hit.localPosition());
546  LocalPoint lp = theChamber->toLocal(gp);
547 
548  // ptc: Local position of hit w.r.t. chamber
549  double u = lp.x();
550  double v = lp.y();
551  double z = lp.z();
552 
553  // ptc: Covariance matrix of local errors
554  CLHEP::HepMatrix IC(2,2);
555  IC(1,1) = hit.localPositionError().xx();
556  IC(1,2) = hit.localPositionError().xy();
557  IC(2,1) = IC(1,2); // since Cov is symmetric
558  IC(2,2) = hit.localPositionError().yy();
559 
560  // ptc: Invert covariance matrix (and trap if it fails!)
561  int ierr = 0;
562  IC.invert(ierr); // inverts in place
563  if (ierr != 0) {
564  LogDebug("CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC << "\n";
565 
566  //@@ NOW WHAT TO DO? Exception? Return? Ignore?
567  //@@ Implement throw
568  }
569 
570  // ptc: Note that IC is no longer Cov but Cov^-1
571  M(1,1) += IC(1,1);
572  M(1,2) += IC(1,2);
573  M(1,3) += IC(1,1) * z;
574  M(1,4) += IC(1,2) * z;
575  B(1) += u * IC(1,1) + v * IC(1,2);
576 
577  M(2,1) += IC(2,1);
578  M(2,2) += IC(2,2);
579  M(2,3) += IC(2,1) * z;
580  M(2,4) += IC(2,2) * z;
581  B(2) += u * IC(2,1) + v * IC(2,2);
582 
583  M(3,1) += IC(1,1) * z;
584  M(3,2) += IC(1,2) * z;
585  M(3,3) += IC(1,1) * z * z;
586  M(3,4) += IC(1,2) * z * z;
587  B(3) += ( u * IC(1,1) + v * IC(1,2) ) * z;
588 
589  M(4,1) += IC(2,1) * z;
590  M(4,2) += IC(2,2) * z;
591  M(4,3) += IC(2,1) * z * z;
592  M(4,4) += IC(2,2) * z * z;
593  B(4) += ( u * IC(2,1) + v * IC(2,2) ) * z;
594  }
595 
596  // Solve the matrix equation using CLHEP's 'solve'
597  //@@ ptc: CAN solve FAIL?? UNCLEAR FROM (LACK OF) CLHEP DOC
598  CLHEP::HepVector p = solve(M, B);
599 
600  // Update member variables uz, vz, theOrigin
601  // Note it has local z = 0
602  theOrigin = LocalPoint(p(1), p(2), 0.);
603  uz = p(3);
604  vz = p(4);
605 }
606 
608 
609  // The chi-squared is (m-Ap)T E (m-Ap)
610  // where T denotes transpose.
611  // This collapses to a simple sum over contributions from each
612  // pair of measurements.
613  float u0 = theOrigin.x();
614  float v0 = theOrigin.y();
615  double chsq = 0.;
616 
617  ChamberHitContainer::const_iterator ih;
618  for (ih = proto_segment.begin(); ih != proto_segment.end(); ++ih) {
619 
620  const CSCRecHit2D& hit = (**ih);
621  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
622  GlobalPoint gp = layer->toGlobal(hit.localPosition());
623  LocalPoint lp = theChamber->toLocal(gp); // FIX !!
624 
625  double hu = lp.x();
626  double hv = lp.y();
627  double hz = lp.z();
628 
629  double du = u0 + uz * hz - hu;
630  double dv = v0 + vz * hz - hv;
631 
632  CLHEP::HepMatrix IC(2,2);
633  IC(1,1) = hit.localPositionError().xx();
634  IC(1,2) = hit.localPositionError().xy();
635  IC(2,1) = IC(1,2);
636  IC(2,2) = hit.localPositionError().yy();
637 
638  // Invert covariance matrix
639  int ierr = 0;
640  IC.invert(ierr);
641  if (ierr != 0) {
642  LogDebug("CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC << "\n";
643 
644  // @@ NOW WHAT TO DO? Exception? Return? Ignore?
645  }
646 
647  chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
648  }
649  theChi2 = chsq;
650 }
651 
653  // Always enforce direction of segment to point from IP outwards
654  // (Incorrect for particles not coming from IP, of course.)
655 
656  double dxdz = uz;
657  double dydz = vz;
658  double dz = 1./sqrt(1. + dxdz*dxdz + dydz*dydz);
659  double dx = dz*dxdz;
660  double dy = dz*dydz;
661  LocalVector localDir(dx,dy,dz);
662 
663  // localDir may need sign flip to ensure it points outward from IP
664  // ptc: Examine its direction and origin in global z: to point outward
665  // the localDir should always have same sign as global z...
666 
667  double globalZpos = ( theChamber->toGlobal( theOrigin ) ).z();
668  double globalZdir = ( theChamber->toGlobal( localDir ) ).z();
669  double directionSign = globalZpos * globalZdir;
670 
671  theDirection = (directionSign * localDir).unit();
672 
673 }
674 
675 bool CSCSegAlgoSK::hasHitOnLayer(int layer) const {
676 
677  // Is there is already a hit on this layer?
679 
680  for(it = proto_segment.begin(); it != proto_segment.end(); it++)
681  if ((*it)->cscDetId().layer() == layer)
682  return true;
683 
684  return false;
685 }
686 
687 bool CSCSegAlgoSK::replaceHit(const CSCRecHit2D* h, int layer) {
688 
689  // replace a hit from a layer
690  ChamberHitContainer::iterator it;
691  for (it = proto_segment.begin(); it != proto_segment.end();) {
692  if ((*it)->cscDetId().layer() == layer)
693  it = proto_segment.erase(it);
694  else
695  ++it;
696  }
697 
698  return addHit(h, layer);
699 }
700 
702 
703  // compare the chi2 of two segments
704  double oldChi2 = theChi2;
705  LocalPoint oldOrigin = theOrigin;
706  LocalVector oldDirection = theDirection;
707  ChamberHitContainer oldSegment = proto_segment;
708 
709  bool ok = replaceHit(h, layer);
710 
711  if (ok) {
712  LogDebug("CSC") << " hit in same layer as a hit on segment; try replacing old one..."
713  << " chi2 new: " << theChi2 << " old: " << oldChi2 << "\n";
714  }
715 
716  if ((theChi2 < oldChi2) && (ok)) {
717  LogDebug("CSC") << " segment with replaced hit is better.\n";
718  }
719  else {
720  proto_segment = oldSegment;
721  theChi2 = oldChi2;
722  theOrigin = oldOrigin;
723  theDirection = oldDirection;
724  }
725 }
726 
728 
729  double oldChi2 = theChi2;
730  LocalPoint oldOrigin = theOrigin;
731  LocalVector oldDirection = theDirection;
732  ChamberHitContainer oldSegment = proto_segment;
733 
734  bool ok = addHit(h, layer);
735 
736  if (ok) {
737  LogDebug("CSC") << " hit in new layer: added to segment, new chi2: "
738  << theChi2 << "\n";
739  }
740 
741  int ndf = 2*proto_segment.size() - 4;
742 
743  if (ok && ((ndf <= 0) || (theChi2/ndf < chi2Max))) {
744  LogDebug("CSC") << " segment with added hit is good.\n" ; // FIX
745  }
746  else {
747  proto_segment = oldSegment;
748  theChi2 = oldChi2;
749  theOrigin = oldOrigin;
750  theDirection = oldDirection;
751  }
752 }
753 
754 CLHEP::HepMatrix CSCSegAlgoSK::derivativeMatrix() const {
755 
756  ChamberHitContainer::const_iterator it;
757  int nhits = proto_segment.size();
758  CLHEP::HepMatrix matrix(2*nhits, 4);
759  int row = 0;
760 
761  for(it = proto_segment.begin(); it != proto_segment.end(); ++it) {
762 
763  const CSCRecHit2D& hit = (**it);
764  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
765  GlobalPoint gp = layer->toGlobal(hit.localPosition());
766  LocalPoint lp = theChamber->toLocal(gp);
767  float z = lp.z();
768  ++row;
769  matrix(row, 1) = 1.;
770  matrix(row, 3) = z;
771  ++row;
772  matrix(row, 2) = 1.;
773  matrix(row, 4) = z;
774  }
775  return matrix;
776 }
777 
778 
780 
781  std::vector<const CSCRecHit2D*>::const_iterator it;
782  int nhits = proto_segment.size();
783  AlgebraicSymMatrix matrix(2*nhits, 0);
784  int row = 0;
785 
786  for (it = proto_segment.begin(); it != proto_segment.end(); ++it) {
787 
788  const CSCRecHit2D& hit = (**it);
789  ++row;
790  matrix(row, row) = hit.localPositionError().xx();
791  matrix(row, row+1) = hit.localPositionError().xy();
792  ++row;
793  matrix(row, row-1) = hit.localPositionError().xy();
794  matrix(row, row) = hit.localPositionError().yy();
795  }
796  int ierr;
797  matrix.invert(ierr);
798  return matrix;
799 }
800 
802 
805 
806  // (AT W A)^-1
807  // from http://www.phys.ufl.edu/~avery/fitting.html, part I
808  int ierr;
809  AlgebraicSymMatrix result = weights.similarityT(A);
810  result.invert(ierr);
811 
812  // blithely assuming the inverting never fails...
813  return result;
814 }
815 
817 
818  // The CSCSegment needs the error matrix re-arranged
819 
820  AlgebraicSymMatrix hold( a );
821 
822  // errors on slopes into upper left
823  a(1,1) = hold(3,3);
824  a(1,2) = hold(3,4);
825  a(2,1) = hold(4,3);
826  a(2,2) = hold(4,4);
827 
828  // errors on positions into lower right
829  a(3,3) = hold(1,1);
830  a(3,4) = hold(1,2);
831  a(4,3) = hold(2,1);
832  a(4,4) = hold(2,2);
833 
834  // off-diagonal elements remain unchanged
835 
836 }
#define LogDebug(id)
void dumpHits(const ChamberHitContainer &rechits) const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
float dPhiFineMax
Definition: CSCSegAlgoSK.h:145
int i
Definition: DBlmapReader.cc:9
float xx() const
Definition: LocalError.h:24
CSCDetId cscDetId() const
Definition: CSCRecHit2D.h:52
AlgebraicSymMatrix weightMatrix(void) const
T perp() const
Definition: PV3DBase.h:71
void fillChiSquared(void)
bool areHitsCloseInLocalX(const CSCRecHit2D *h1, const CSCRecHit2D *h2) const
Utility functions.
ChamberHitContainer proto_segment
Definition: CSCSegAlgoSK.h:134
void flagHitsAsUsed(const ChamberHitContainer &rechitsInChamber, BoolContainer &used) const
std::vector< const CSCRecHit2D * >::const_iterator ChamberHitContainerCIt
Definition: CSCSegAlgoSK.h:53
double_binary B
Definition: DDStreamer.cc:234
bool replaceHit(const CSCRecHit2D *h, int layer)
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
Geom::Phi< T > phi() const
Definition: PV3DBase.h:68
T y() const
Definition: PV3DBase.h:62
#define abs(x)
Definition: mlp_lapack.h:159
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:62
LocalPoint theOrigin
Definition: CSCSegAlgoSK.h:138
std::vector< CSCSegment > buildSegments(ChamberHitContainer rechits)
Definition: CSCSegAlgoSK.cc:52
int layer() const
Definition: CSCDetId.h:63
bool addHit(const CSCRecHit2D *hit, int layer)
Utility functions.
double double double z
void fitSlopes(void)
bool isHitNearSegment(const CSCRecHit2D *h) const
std::string chamberTypeName() const
float windowScale
Definition: CSCSegAlgoSK.h:141
AlgebraicSymMatrix calculateError(void) const
float xy() const
Definition: LocalError.h:25
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:41
const CSCChamberSpecs * specs() const
Definition: CSCChamber.h:42
CLHEP::HepMatrix AlgebraicMatrix
LocalError localPositionError() const
Definition: CSCRecHit2D.h:51
float yy() const
Definition: LocalError.h:26
void increaseProtoSegment(const CSCRecHit2D *h, int layer)
string unit
Definition: csvLumiCalc.py:46
T sqrt(T t)
Definition: SSEVec.h:46
T z() const
Definition: PV3DBase.h:63
tuple result
Definition: query.py:137
void flipErrors(AlgebraicSymMatrix &) const
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:41
double f[11][100]
float phiAtZ(float z) const
void updateParameters(void)
LocalVector theDirection
Definition: CSCSegAlgoSK.h:139
CSCSegAlgoSK(const edm::ParameterSet &ps)
Constructor.
Definition: CSCSegAlgoSK.cc:23
#define M_PI
Definition: BFit3D.cc:3
bool hasHitOnLayer(int layer) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
void tryAddingHitsToSegment(const ChamberHitContainer &rechitsInChamber, BoolContainer used, LayerIndex layerIndex, const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2)
const std::string myName
Definition: CSCSegAlgoSK.h:135
std::vector< int > LayerIndex
Typedefs.
Definition: CSCSegAlgoSK.h:51
float dRPhiFineMax
Definition: CSCSegAlgoSK.h:144
double a
Definition: hdecay.h:121
CLHEP::HepSymMatrix AlgebraicSymMatrix
bool isSegmentGood(const ChamberHitContainer &rechitsInChamber) const
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
double theChi2
Definition: CSCSegAlgoSK.h:137
std::vector< const CSCRecHit2D * > ChamberHitContainer
Definition: CSCSegAlgoSK.h:52
std::vector< CSCSegment > run(const CSCChamber *aChamber, ChamberHitContainer rechits)
Definition: CSCSegAlgoSK.cc:47
CLHEP::HepMatrix derivativeMatrix(void) const
Definition: DDAxes.h:10
LocalPoint localPosition() const
Definition: CSCRecHit2D.h:50
T x() const
Definition: PV3DBase.h:61
bool areHitsCloseInGlobalPhi(const CSCRecHit2D *h1, const CSCRecHit2D *h2) const
mathSSE::Vec4< T > v
std::deque< bool > BoolContainer
Definition: CSCSegAlgoSK.h:60
void compareProtoSegment(const CSCRecHit2D *h, int layer)
const CSCChamber * theChamber
Definition: CSCSegAlgoSK.h:133
Definition: DDAxes.h:10
void fillLocalDirection(void)