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