CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoLocalMuon/CSCSegment/src/CSCSegAlgoSK.cc

Go to the documentation of this file.
00001 
00009 #include "CSCSegAlgoSK.h"
00010 #include "Geometry/CSCGeometry/interface/CSCLayer.h"
00011 // For clhep Matrix::solve
00012 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
00013 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00014 
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 
00018 #include <algorithm>
00019 #include <cmath>
00020 #include <iostream>
00021 #include <string>
00022 
00023 CSCSegAlgoSK::CSCSegAlgoSK(const edm::ParameterSet& ps) : CSCSegmentAlgorithm(ps),
00024         myName("CSCSegAlgoSK") {
00025         
00026   debugInfo = ps.getUntrackedParameter<bool>("verboseInfo");
00027     
00028   dRPhiMax       = ps.getParameter<double>("dRPhiMax");
00029   dPhiMax        = ps.getParameter<double>("dPhiMax");
00030   dRPhiFineMax   = ps.getParameter<double>("dRPhiFineMax");
00031   dPhiFineMax    = ps.getParameter<double>("dPhiFineMax");
00032   chi2Max        = ps.getParameter<double>("chi2Max");
00033   wideSeg        = ps.getParameter<double>("wideSeg");
00034   minLayersApart = ps.getParameter<int>("minLayersApart");
00035         
00036   LogDebug("CSC") << myName << " has algorithm cuts set to: \n"
00037                   << "--------------------------------------------------------------------\n"
00038                   << "dRPhiMax     = " << dRPhiMax << '\n'
00039                   << "dPhiMax      = " << dPhiMax << '\n'
00040                   << "dRPhiFineMax = " << dRPhiFineMax << '\n'
00041                   << "dPhiFineMax  = " << dPhiFineMax << '\n'
00042                   << "chi2Max      = " << chi2Max << '\n'
00043                   << "wideSeg      = " << wideSeg << '\n'
00044                   << "minLayersApart = " << minLayersApart << std::endl;
00045 }
00046 
00047 std::vector<CSCSegment> CSCSegAlgoSK::run(const CSCChamber* aChamber, ChamberHitContainer rechits) {
00048     theChamber = aChamber; 
00049     return buildSegments(rechits); 
00050 }
00051 
00052 std::vector<CSCSegment> CSCSegAlgoSK::buildSegments(ChamberHitContainer rechits) {
00053         
00054   LogDebug("CSC") << "*********************************************";
00055   LogDebug("CSC") << "Start segment building in the new chamber: " << theChamber->specs()->chamberTypeName();
00056   LogDebug("CSC") << "*********************************************";
00057   
00058   
00059   LayerIndex layerIndex(rechits.size());
00060   
00061   for(unsigned int i = 0; i < rechits.size(); i++) {
00062     
00063     layerIndex[i] = rechits[i]->cscDetId().layer();
00064   }
00065   
00066   double z1 = theChamber->layer(1)->position().z();
00067   double z6 = theChamber->layer(6)->position().z();
00068   
00069   if ( z1 > 0. ) {
00070     if ( z1 > z6 ) { 
00071       reverse(layerIndex.begin(), layerIndex.end());
00072       reverse(rechits.begin(), rechits.end());
00073     }    
00074   }
00075   else if ( z1 < 0. ) {
00076     if ( z1 < z6 ) {
00077       reverse(layerIndex.begin(), layerIndex.end());
00078       reverse(rechits.begin(), rechits.end());
00079     }    
00080   }
00081   
00082   if (debugInfo) {
00083     // dump after sorting
00084     dumpHits(rechits);
00085   }
00086   
00087   if (rechits.size() < 2) {
00088     LogDebug("CSC") << myName << ": " << rechits.size() << 
00089       "  hit(s) in chamber is not enough to build a segment.\n";
00090     return std::vector<CSCSegment>(); 
00091   }
00092         
00093   // We have at least 2 hits. We intend to try all possible pairs of hits to start 
00094   // segment building. 'All possible' means each hit lies on different layers in the chamber.
00095   // BUT... once a hit has been assigned to a segment, we don't consider
00096   // it again.
00097   
00098   // Choose first hit (as close to IP as possible) h1 and 
00099   // second hit (as far from IP as possible) h2
00100   // To do this we iterate over hits in the chamber by layer - pick two layers.
00101   // @@ Require the two layers are at least 3 layers apart. May need tuning?
00102   // Then we iterate over hits within each of these layers and pick h1 and h2 from these.
00103   // If they are 'close enough' we build an empty segment.
00104   // Then try adding hits to this segment.
00105   
00106   // Initialize flags that a given hit has been allocated to a segment
00107   BoolContainer used(rechits.size(), false);
00108   
00109   // Define buffer for segments we build 
00110   std::vector<CSCSegment> segments;
00111   
00112   ChamberHitContainerCIt ib = rechits.begin();
00113   ChamberHitContainerCIt ie = rechits.end();
00114   
00115   // Possibly allow 2 passes, second widening scale factor for cuts
00116   windowScale = 1.; // scale factor for cuts
00117   
00118   int npass = (wideSeg > 1.)? 2 : 1;
00119   
00120   for (int ipass = 0; ipass < npass; ++ipass) {
00121     for (ChamberHitContainerCIt i1 = ib; i1 != ie; ++i1) {
00122       bool segok = false;
00123       if(used[i1-ib]) 
00124         continue;
00125       
00126       int layer1 = layerIndex[i1-ib]; //(*i1)->cscDetId().layer();
00127       const CSCRecHit2D* h1 = *i1;
00128       
00129       for (ChamberHitContainerCIt i2 = ie-1; i2 != i1; --i2) {
00130         if(used[i2-ib]) 
00131           continue;
00132         
00133         int layer2 = layerIndex[i2-ib]; //(*i2)->cscDetId().layer();
00134                                 
00135         if (abs(layer2 - layer1) < minLayersApart) 
00136           break;
00137         const CSCRecHit2D* h2 = *i2;
00138         
00139         if (areHitsCloseInLocalX(h1, h2) && areHitsCloseInGlobalPhi(h1, h2)) {
00140           
00141           proto_segment.clear();
00142           
00143           const CSCLayer* l1 = theChamber->layer(layer1);
00144           GlobalPoint gp1 = l1->toGlobal(h1->localPosition());                                  
00145           const CSCLayer* l2 = theChamber->layer(layer2);
00146           GlobalPoint gp2 = l2->toGlobal(h2->localPosition());                                  
00147           LogDebug("CSC") << "start new segment from hits " << "h1: " 
00148                           << gp1 << " - h2: " << gp2 << "\n";
00149           
00150           if (!addHit(h1, layer1)) { 
00151             LogDebug("CSC") << "  failed to add hit h1\n";
00152             continue;
00153           }
00154           
00155           if (!addHit(h2, layer2)) { 
00156             LogDebug("CSC") << "  failed to add hit h2\n";
00157             continue;
00158           }
00159           
00160           tryAddingHitsToSegment(rechits, used, layerIndex, i1, i2); 
00161           
00162           // Check no. of hits on segment, and if enough flag them as used
00163           // and store the segment
00164           segok = isSegmentGood(rechits);
00165           if (segok) {
00166             flagHitsAsUsed(rechits, used);
00167             // Copy the proto_segment and its properties inside a CSCSegment.
00168             // Then fill the segment vector..
00169             
00170             if (proto_segment.empty()) {
00171               LogDebug("CSC") << "No segment has been found !!!\n";
00172             }   
00173             else {
00174               
00175               // calculate error matrix...
00176               AlgebraicSymMatrix errors = calculateError();
00177 
00178               // but reorder components to match what's required by TrackingRecHit interface
00179               // i.e. slopes first, then positions
00180 
00181               flipErrors( errors );
00182 
00183               CSCSegment temp(proto_segment, theOrigin, theDirection, errors, theChi2); 
00184               
00185               LogDebug("CSC") << "Found a segment !!!\n";
00186               segments.push_back(temp); 
00187             }
00188           }
00189         }  //   h1 & h2 close
00190         
00191         if (segok) 
00192           break;
00193       }  //  i2
00194     }  //  i1
00195     
00196     if (segments.size() > 1)  
00197       break;  // only change window if no segments found
00198     
00199     // Increase cut windows by factor of wideSeg
00200     windowScale = wideSeg;
00201     
00202   }  //  ipass
00203   
00204   // Give the segments to the CSCChamber
00205   return segments;
00206 }
00207 
00208 void CSCSegAlgoSK::tryAddingHitsToSegment(const ChamberHitContainer& rechits, 
00209                                           BoolContainer used, LayerIndex layerIndex,
00210                                           const ChamberHitContainerCIt i1, const ChamberHitContainerCIt i2) {
00211 
00212   // Iterate over the layers with hits in the chamber
00213   // Skip the layers containing the segment endpoints
00214   // Test each hit on the other layers to see if it is near the segment
00215   // If it is, see whether there is already a hit on the segment from the same layer
00216   //    - if so, and there are more than 2 hits on the segment, copy the segment,
00217   //      replace the old hit with the new hit. If the new segment chi2 is better
00218   //      then replace the original segment with the new one (by swap)
00219   //    - if not, copy the segment, add the hit. If the new chi2/dof is still satisfactory
00220   //      then replace the original segment with the new one (by swap)
00221   
00222   ChamberHitContainerCIt ib = rechits.begin();
00223   ChamberHitContainerCIt ie = rechits.end();
00224   
00225   for (ChamberHitContainerCIt i = ib; i != ie; ++i) {
00226     if (i == i1 || i == i2 || used[i-ib]) 
00227       continue; 
00228     
00229     int layer = layerIndex[i-ib]; 
00230     const CSCRecHit2D* h = *i;
00231     if (isHitNearSegment(h)) {
00232       
00233       GlobalPoint gp1 = theChamber->layer(layer)->toGlobal(h->localPosition());         
00234       LogDebug("CSC") << "    hit at global " << gp1 << " is near segment\n.";
00235       
00236       // Don't consider alternate hits on layers holding the two starting points
00237       if (hasHitOnLayer(layer)) {
00238         if (proto_segment.size() <= 2) {
00239           LogDebug("CSC") << "    " << proto_segment.size() 
00240                           << " hits on segment...skip hit on same layer.\n";
00241           continue;
00242         }
00243         compareProtoSegment(h, layer);
00244       } 
00245       else                      
00246         increaseProtoSegment(h, layer);
00247     }   // h & seg close
00248   }   // i
00249 }
00250 
00251 bool CSCSegAlgoSK::areHitsCloseInLocalX(const CSCRecHit2D* h1, const CSCRecHit2D* h2) const {
00252   float h1x = h1->localPosition().x();
00253   float h2x = h2->localPosition().x();
00254   float deltaX = (h1->localPosition()-h2->localPosition()).x();
00255   LogDebug("CSC") << "    Hits at local x= " << h1x << ", " 
00256                   << h2x << " have separation= " << deltaX;
00257   return (fabs(deltaX) < (dRPhiMax * windowScale))? true:false;   // +v
00258 }
00259 
00260 bool CSCSegAlgoSK::areHitsCloseInGlobalPhi(const CSCRecHit2D* h1, const CSCRecHit2D* h2) const {
00261   
00262   const CSCLayer* l1 = theChamber->layer(h1->cscDetId().layer());
00263   GlobalPoint gp1 = l1->toGlobal(h1->localPosition());  
00264   const CSCLayer* l2 = theChamber->layer(h2->cscDetId().layer());
00265   GlobalPoint gp2 = l2->toGlobal(h2->localPosition());  
00266   
00267   float h1p = gp1.phi();
00268   float h2p = gp2.phi();
00269   float dphi12 = h1p - h2p;
00270   
00271   // Into range [-pi, pi) (phi() returns values in this range)
00272   if (dphi12 < -M_PI) 
00273     dphi12 += 2.*M_PI;  
00274   if (dphi12 >  M_PI) 
00275     dphi12 -= 2.*M_PI;
00276   LogDebug("CSC") << "    Hits at global phi= " << h1p << ", " 
00277                   << h2p << " have separation= " << dphi12;
00278   return (fabs(dphi12) < (dPhiMax * windowScale))? true:false;  // +v
00279 }
00280 
00281 bool CSCSegAlgoSK::isHitNearSegment(const CSCRecHit2D* h) const {
00282   
00283   // Is hit near segment? 
00284   // Requires deltaphi and rxy*deltaphi within ranges specified
00285   // in parameter set, where rxy=sqrt(x**2+y**2) of hit itself.
00286   // Note that to make intuitive cuts on delta(phi) one must work in
00287   // phi range (-pi, +pi] not [0, 2pi)
00288   
00289   const CSCLayer* l1 = theChamber->layer(h->cscDetId().layer());
00290   GlobalPoint hp = l1->toGlobal(h->localPosition());    
00291   
00292   float hphi = hp.phi();          // in (-pi, +pi]
00293   if (hphi < 0.) 
00294     hphi += 2.*M_PI;            // into range [0, 2pi)
00295   float sphi = phiAtZ(hp.z());    // in [0, 2*pi)
00296   float phidif = sphi-hphi;
00297   if (phidif < 0.) 
00298     phidif += 2.*M_PI;          // into range [0, 2pi)
00299   if (phidif > M_PI) 
00300     phidif -= 2.*M_PI;          // into range (-pi, pi]
00301   
00302   float dRPhi = fabs(phidif)*hp.perp();
00303   LogDebug("CSC") << "    is hit at phi_h= " << hphi << " near segment phi_seg= " << sphi 
00304                   << "? is " << dRPhi << "<" << dRPhiFineMax*windowScale << " ? " 
00305                   << " and is |" << phidif << "|<" << dPhiFineMax*windowScale << " ?";
00306   
00307   return ((dRPhi < dRPhiFineMax*windowScale) && 
00308           (fabs(phidif) < dPhiFineMax*windowScale))? true:false;  // +v
00309 }
00310 
00311 float CSCSegAlgoSK::phiAtZ(float z) const {
00312   
00313   // Returns a phi in [ 0, 2*pi )
00314   const CSCLayer* l1 = theChamber->layer((*(proto_segment.begin()))->cscDetId().layer());
00315   GlobalPoint gp = l1->toGlobal(theOrigin);     
00316   GlobalVector gv = l1->toGlobal(theDirection); 
00317   
00318   float x = gp.x() + (gv.x()/gv.z())*(z - gp.z());
00319   float y = gp.y() + (gv.y()/gv.z())*(z - gp.z());
00320   float phi = atan2(y, x);
00321   if (phi < 0.f ) 
00322     phi += 2. * M_PI;
00323   
00324   return phi ;
00325 }
00326 
00327 void CSCSegAlgoSK::dumpHits(const ChamberHitContainer& rechits) const {
00328   
00329   // Dump positions of RecHit's in each CSCChamber
00330   ChamberHitContainerCIt it;
00331   edm::LogInfo("CSC") << "CSCChamber rechit dump.\n";   
00332   for(it=rechits.begin(); it!=rechits.end(); it++) {
00333     
00334     const CSCLayer* l1 = theChamber->layer((*it)->cscDetId().layer());
00335     GlobalPoint gp1 = l1->toGlobal((*it)->localPosition());     
00336     
00337     edm::LogInfo("CSC") << "Global pos.: " << gp1 << ", phi: " << gp1.phi() << ". Local position: "
00338                         << (*it)->localPosition() << ", phi: "
00339                         << (*it)->localPosition().phi() << ". Layer: "
00340                         << (*it)->cscDetId().layer() << "\n";
00341   }     
00342 }
00343 
00344 bool CSCSegAlgoSK::isSegmentGood(const ChamberHitContainer& rechitsInChamber) const {
00345   
00346   // If the chamber has 20 hits or fewer, require at least 3 hits on segment
00347   // If the chamber has >20 hits require at least 4 hits
00348   //@@ THESE VALUES SHOULD BECOME PARAMETERS?
00349   bool ok = false;
00350   
00351   unsigned int iadd = ( rechitsInChamber.size() > 20 )?  1 : 0;  
00352   
00353   if (windowScale > 1.)
00354     iadd = 1;
00355   
00356   if (proto_segment.size() >= 3+iadd)
00357     ok = true;
00358   
00359   return ok;
00360 }
00361 
00362 void CSCSegAlgoSK::flagHitsAsUsed(const ChamberHitContainer& rechitsInChamber, 
00363                                   BoolContainer& used ) const {
00364   
00365   // Flag hits on segment as used
00366   ChamberHitContainerCIt ib = rechitsInChamber.begin();
00367   ChamberHitContainerCIt hi, iu;
00368   
00369   for(hi = proto_segment.begin(); hi != proto_segment.end(); ++hi) {
00370     for(iu = ib; iu != rechitsInChamber.end(); ++iu) {
00371       if(*hi == *iu)
00372         used[iu-ib] = true;
00373     }
00374   }
00375 }
00376 
00377 bool CSCSegAlgoSK::addHit(const CSCRecHit2D* aHit, int layer) {
00378   
00379   // Return true if hit was added successfully 
00380   // (and then parameters are updated).
00381   // Return false if there is already a hit on the same layer, or insert failed.
00382   
00383   ChamberHitContainer::const_iterator it;
00384   
00385   for(it = proto_segment.begin(); it != proto_segment.end(); it++)
00386     if (((*it)->cscDetId().layer() == layer) && (aHit != (*it)))
00387       return false; 
00388   
00389   proto_segment.push_back(aHit);
00390   updateParameters();
00391   
00392   return true;
00393 }
00394 
00395 void CSCSegAlgoSK::updateParameters() {
00396   
00397   // Note that we need local position of a RecHit w.r.t. the CHAMBER
00398   // and the RecHit itself only knows its local position w.r.t.
00399   // the LAYER, so need to explicitly transform to global position.
00400   
00401   //  no. of hits in the RecHitsOnSegment
00402   //  By construction this is the no. of layers with hits
00403   //  since we allow just one hit per layer in a segment.
00404   
00405   int nh = proto_segment.size();
00406   
00407   // First hit added to a segment must always fail here
00408   if (nh < 2) 
00409     return;
00410   
00411   if (nh == 2) {
00412     
00413     // Once we have two hits we can calculate a straight line 
00414     // (or rather, a straight line for each projection in xz and yz.)
00415     ChamberHitContainer::const_iterator ih = proto_segment.begin();
00416     int il1 = (*ih)->cscDetId().layer();
00417     const CSCRecHit2D& h1 = (**ih);
00418     ++ih;    
00419     int il2 = (*ih)->cscDetId().layer();
00420     const CSCRecHit2D& h2 = (**ih);
00421     
00422     //@@ Skip if on same layer, but should be impossible
00423     if (il1 == il2) 
00424       return;
00425     
00426     const CSCLayer* layer1 = theChamber->layer(il1);
00427     const CSCLayer* layer2 = theChamber->layer(il2);
00428     
00429     GlobalPoint h1glopos = layer1->toGlobal(h1.localPosition());
00430     GlobalPoint h2glopos = layer2->toGlobal(h2.localPosition());
00431     
00432     // localPosition is position of hit wrt layer (so local z = 0)
00433     theOrigin = h1.localPosition();
00434     
00435     // We want hit wrt chamber (and local z will be != 0)
00436     LocalPoint h1pos = theChamber->toLocal(h1glopos);  
00437     LocalPoint h2pos = theChamber->toLocal(h2glopos);  
00438     
00439     float dz = h2pos.z()-h1pos.z();
00440     uz = (h2pos.x()-h1pos.x())/dz ;
00441     vz = (h2pos.y()-h1pos.y())/dz ;
00442     
00443     theChi2 = 0.;
00444   }
00445   else if (nh > 2) {
00446     
00447     // When we have more than two hits then we can fit projections to straight lines
00448     fitSlopes();  
00449     fillChiSquared();
00450   } // end of 'if' testing no. of hits
00451   
00452   fillLocalDirection(); 
00453 }
00454 
00455 void CSCSegAlgoSK::fitSlopes() {
00456   
00457   // Update parameters of fit
00458   // ptc 13-Aug-02: This does a linear least-squares fit
00459   // to the hits associated with the segment, in the z projection.
00460   
00461   // In principle perhaps one could fit the strip and wire
00462   // measurements (u, v respectively), to
00463   // u = u0 + uz * z
00464   // v = v0 + vz * z
00465   // where u0, uz, v0, vz are the parameters resulting from the fit.
00466   // But what is actually done is fit to the local x, y coordinates 
00467   // of the RecHits. However the strip measurement controls the precision
00468   // of x, and the wire measurement controls that of y.
00469   // Precision in local coordinate:
00470   //       u (strip, sigma~200um), v (wire, sigma~1cm)
00471   
00472   // I have verified that this code agrees with the formulation given
00473   // on p246-247 of 'Data analysis techniques for high-energy physics
00474   // experiments' by Bock, Grote, Notz & Regler, and that on p111-113
00475   // of 'Statistics' by Barlow.
00476   
00477   // Formulate the matrix equation representing the least-squares fit
00478   // We have a vector of measurements m, which is a 2n x 1 dim matrix
00479   // The transpose mT is (u1, v1, u2, v2, ..., un, vn)
00480   // where ui is the strip-associated measurement and vi is the
00481   // wire-associated measurement for a given RecHit i.
00482   // The fit is to
00483   // u = u0 + uz * z
00484   // v = v0 + vz * z
00485   // where u0, uz, v0, vz are the parameters to be obtained from the fit.
00486   // These are contained in a vector p which is a 4x1 dim matrix, and
00487   // its transpose pT is (u0, v0, uz, vz). Note the ordering!
00488   // The covariance matrix for each pair of measurements is 2 x 2 and
00489   // the inverse of this is the error matrix E.
00490   // The error matrix for the whole set of n measurements is a diagonal
00491   // matrix with diagonal elements the individual 2 x 2 error matrices
00492   // (because the inverse of a diagonal matrix is a diagonal matrix
00493   // with each element the inverse of the original.)
00494   
00495   // The matrix 'matrix' in method 'CSCSegment::weightMatrix()' is this 
00496   // block-diagonal overall covariance matrix. It is inverted to the 
00497   // block-diagonal error matrix right before it is returned.
00498   
00499   // Use the matrix A defined by
00500   //    1   0   z1  0
00501   //    0   1   0   z1
00502   //    1   0   z2  0
00503   //    0   1   0   z2
00504   //    ..  ..  ..  ..
00505   //    1   0   zn  0
00506   //    0   1   0   zn
00507   
00508   // The matrix A is returned by 'CSCSegment::derivativeMatrix()'.
00509   
00510   // Then the normal equations are encapsulated in the matrix equation
00511   //
00512   //    (AT E A)p = (AT E)m
00513   //
00514   // where AT is the transpose of A.
00515   // We'll call the combined matrix on the LHS, M, and that on the RHS, B:
00516   //     M p = B
00517   
00518   // We solve this for the parameter vector, p.
00519   // The elements of M and B then involve sums over the hits
00520   
00521   // The error matrix of the parameters is obtained by 
00522   // (AT E A)^-1 calculated in 'calculateError()'.
00523   
00524   // The 4 values in p can be accessed from 'CSCSegment::parameters()'
00525   // in the order uz, vz, u0, v0.
00526   
00527   // NOTE 1
00528   // Do the #hits = 2 case separately.
00529   // (I hope they're not on the same layer! They should not be, by construction.)
00530   
00531   // NOTE 2
00532   // We need local position of a RecHit w.r.t. the CHAMBER
00533   // and the RecHit itself only knows its local position w.r.t.
00534   // the LAYER, so we must explicitly transform global position.
00535   
00536   CLHEP::HepMatrix M(4,4,0);
00537   CLHEP::HepVector B(4,0);
00538   
00539   ChamberHitContainer::const_iterator ih = proto_segment.begin();
00540   
00541   for (ih = proto_segment.begin(); ih != proto_segment.end(); ++ih) {
00542     
00543     const CSCRecHit2D& hit = (**ih);
00544     const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
00545     GlobalPoint gp = layer->toGlobal(hit.localPosition());
00546     LocalPoint  lp  = theChamber->toLocal(gp); 
00547     
00548     // ptc: Local position of hit w.r.t. chamber
00549     double u = lp.x();
00550     double v = lp.y();
00551     double z = lp.z();
00552     
00553     // ptc: Covariance matrix of local errors 
00554     CLHEP::HepMatrix IC(2,2);
00555     IC(1,1) = hit.localPositionError().xx();
00556     IC(1,2) = hit.localPositionError().xy();
00557     IC(2,1) = IC(1,2); // since Cov is symmetric
00558     IC(2,2) = hit.localPositionError().yy();
00559     
00560     // ptc: Invert covariance matrix (and trap if it fails!)
00561     int ierr = 0;
00562     IC.invert(ierr); // inverts in place
00563     if (ierr != 0) {
00564       LogDebug("CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC << "\n";
00565       
00566       //@@ NOW WHAT TO DO? Exception? Return? Ignore?
00567       //@@ Implement throw
00568     }
00569     
00570     // ptc: Note that IC is no longer Cov but Cov^-1
00571     M(1,1) += IC(1,1);
00572     M(1,2) += IC(1,2);
00573     M(1,3) += IC(1,1) * z;
00574     M(1,4) += IC(1,2) * z;
00575     B(1) += u * IC(1,1) + v * IC(1,2);
00576     
00577     M(2,1) += IC(2,1);
00578     M(2,2) += IC(2,2);
00579     M(2,3) += IC(2,1) * z;
00580     M(2,4) += IC(2,2) * z;
00581     B(2) += u * IC(2,1) + v * IC(2,2);
00582     
00583     M(3,1) += IC(1,1) * z;
00584     M(3,2) += IC(1,2) * z;
00585     M(3,3) += IC(1,1) * z * z;
00586     M(3,4) += IC(1,2) * z * z;
00587     B(3) += ( u * IC(1,1) + v * IC(1,2) ) * z;
00588     
00589     M(4,1) += IC(2,1) * z;
00590     M(4,2) += IC(2,2) * z;
00591     M(4,3) += IC(2,1) * z * z;
00592     M(4,4) += IC(2,2) * z * z;
00593     B(4) += ( u * IC(2,1) + v * IC(2,2) ) * z;
00594   }
00595   
00596   // Solve the matrix equation using CLHEP's 'solve'
00597   //@@ ptc: CAN solve FAIL?? UNCLEAR FROM (LACK OF) CLHEP DOC
00598   CLHEP::HepVector p = solve(M, B);
00599   
00600   // Update member variables uz, vz, theOrigin
00601   // Note it has local z = 0
00602   theOrigin = LocalPoint(p(1), p(2), 0.);
00603   uz = p(3);
00604   vz = p(4);
00605 }
00606 
00607 void CSCSegAlgoSK::fillChiSquared() {
00608   
00609   // The chi-squared is (m-Ap)T E (m-Ap)
00610   // where T denotes transpose.
00611   // This collapses to a simple sum over contributions from each
00612   // pair of measurements.
00613   float u0 = theOrigin.x();
00614   float v0 = theOrigin.y();
00615   double chsq = 0.;
00616   
00617   ChamberHitContainer::const_iterator ih;
00618   for (ih = proto_segment.begin(); ih != proto_segment.end(); ++ih) {
00619     
00620     const CSCRecHit2D& hit = (**ih);
00621     const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
00622     GlobalPoint gp = layer->toGlobal(hit.localPosition());
00623     LocalPoint lp = theChamber->toLocal(gp);  // FIX !!
00624     
00625     double hu = lp.x();
00626     double hv = lp.y();
00627     double hz = lp.z();
00628     
00629     double du = u0 + uz * hz - hu;
00630     double dv = v0 + vz * hz - hv;
00631     
00632     CLHEP::HepMatrix IC(2,2);
00633     IC(1,1) = hit.localPositionError().xx();
00634     IC(1,2) = hit.localPositionError().xy();
00635     IC(2,1) = IC(1,2);
00636     IC(2,2) = hit.localPositionError().yy();
00637     
00638     // Invert covariance matrix
00639     int ierr = 0;
00640     IC.invert(ierr);
00641     if (ierr != 0) {
00642       LogDebug("CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC << "\n";
00643       
00644       // @@ NOW WHAT TO DO? Exception? Return? Ignore?
00645     }
00646     
00647     chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
00648   }
00649   theChi2 = chsq;
00650 }
00651 
00652 void CSCSegAlgoSK::fillLocalDirection() {
00653   // Always enforce direction of segment to point from IP outwards
00654   // (Incorrect for particles not coming from IP, of course.)
00655   
00656   double dxdz = uz;
00657   double dydz = vz;
00658   double dz = 1./sqrt(1. + dxdz*dxdz + dydz*dydz);
00659   double dx = dz*dxdz;
00660   double dy = dz*dydz;
00661   LocalVector localDir(dx,dy,dz);
00662   
00663   // localDir may need sign flip to ensure it points outward from IP
00664   // ptc: Examine its direction and origin in global z: to point outward
00665   // the localDir should always have same sign as global z...
00666   
00667   double globalZpos = ( theChamber->toGlobal( theOrigin ) ).z();
00668   double globalZdir = ( theChamber->toGlobal( localDir ) ).z();
00669   double directionSign = globalZpos * globalZdir;
00670   
00671   theDirection = (directionSign * localDir).unit();
00672   
00673 }
00674 
00675 bool CSCSegAlgoSK::hasHitOnLayer(int layer) const {
00676   
00677   // Is there is already a hit on this layer?
00678   ChamberHitContainerCIt it;
00679   
00680   for(it = proto_segment.begin(); it != proto_segment.end(); it++)
00681     if ((*it)->cscDetId().layer() == layer)
00682       return true; 
00683   
00684   return false;
00685 }
00686 
00687 bool CSCSegAlgoSK::replaceHit(const CSCRecHit2D* h, int layer) {
00688   
00689   // replace a hit from a layer 
00690   ChamberHitContainer::iterator it;
00691   for (it = proto_segment.begin(); it != proto_segment.end();) {
00692     if ((*it)->cscDetId().layer() == layer)
00693       it = proto_segment.erase(it);
00694     else
00695       ++it;   
00696   }
00697   
00698   return addHit(h, layer);                                  
00699 }
00700 
00701 void CSCSegAlgoSK::compareProtoSegment(const CSCRecHit2D* h, int layer) {
00702   
00703   // compare the chi2 of two segments
00704   double oldChi2 = theChi2;
00705   LocalPoint oldOrigin = theOrigin;
00706   LocalVector oldDirection = theDirection;
00707   ChamberHitContainer oldSegment = proto_segment;
00708   
00709   bool ok = replaceHit(h, layer);
00710   
00711   if (ok) {
00712     LogDebug("CSC") << "    hit in same layer as a hit on segment; try replacing old one..." 
00713                     << " chi2 new: " << theChi2 << " old: " << oldChi2 << "\n";
00714   }
00715   
00716   if ((theChi2 < oldChi2) && (ok)) {
00717     LogDebug("CSC")  << "    segment with replaced hit is better.\n";
00718   }
00719   else {
00720     proto_segment = oldSegment;
00721     theChi2 = oldChi2;
00722     theOrigin = oldOrigin;
00723     theDirection = oldDirection;
00724   }
00725 }
00726 
00727 void CSCSegAlgoSK::increaseProtoSegment(const CSCRecHit2D* h, int layer) {
00728   
00729   double oldChi2 = theChi2;
00730   LocalPoint oldOrigin = theOrigin;
00731   LocalVector oldDirection = theDirection;
00732   ChamberHitContainer oldSegment = proto_segment;
00733   
00734   bool ok = addHit(h, layer);
00735   
00736   if (ok) {
00737     LogDebug("CSC") << "    hit in new layer: added to segment, new chi2: " 
00738                     << theChi2 << "\n";
00739   }
00740   
00741   int ndf = 2*proto_segment.size() - 4;
00742   
00743   if (ok && ((ndf <= 0) || (theChi2/ndf < chi2Max))) {
00744     LogDebug("CSC") << "    segment with added hit is good.\n" ;                // FIX
00745   }     
00746   else {
00747     proto_segment = oldSegment;
00748     theChi2 = oldChi2;
00749     theOrigin = oldOrigin;
00750     theDirection = oldDirection;
00751   }                     
00752 }               
00753 
00754 CLHEP::HepMatrix CSCSegAlgoSK::derivativeMatrix() const {
00755   
00756   ChamberHitContainer::const_iterator it;
00757   int nhits = proto_segment.size();
00758   CLHEP::HepMatrix matrix(2*nhits, 4);
00759   int row = 0;
00760   
00761   for(it = proto_segment.begin(); it != proto_segment.end(); ++it) {
00762     
00763     const CSCRecHit2D& hit = (**it);
00764     const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
00765     GlobalPoint gp = layer->toGlobal(hit.localPosition());      
00766     LocalPoint lp = theChamber->toLocal(gp); 
00767     float z = lp.z();
00768     ++row;
00769     matrix(row, 1) = 1.;
00770     matrix(row, 3) = z;
00771     ++row;
00772     matrix(row, 2) = 1.;
00773     matrix(row, 4) = z;
00774   }
00775   return matrix;
00776 }
00777 
00778 
00779 AlgebraicSymMatrix CSCSegAlgoSK::weightMatrix() const {
00780   
00781   std::vector<const CSCRecHit2D*>::const_iterator it;
00782   int nhits = proto_segment.size();
00783   AlgebraicSymMatrix matrix(2*nhits, 0);
00784   int row = 0;
00785   
00786   for (it = proto_segment.begin(); it != proto_segment.end(); ++it) {
00787     
00788     const CSCRecHit2D& hit = (**it);
00789     ++row;
00790     matrix(row, row)   = hit.localPositionError().xx();
00791     matrix(row, row+1) = hit.localPositionError().xy();
00792     ++row;
00793     matrix(row, row-1) = hit.localPositionError().xy();
00794     matrix(row, row)   = hit.localPositionError().yy();
00795   }
00796   int ierr;
00797   matrix.invert(ierr);
00798   return matrix;
00799 }
00800 
00801 AlgebraicSymMatrix CSCSegAlgoSK::calculateError() const {
00802   
00803   AlgebraicSymMatrix weights = weightMatrix();
00804   AlgebraicMatrix A = derivativeMatrix();
00805   
00806   // (AT W A)^-1
00807   // from https://www.phys.ufl.edu/~avery/fitting.html, part I
00808   int ierr;
00809   AlgebraicSymMatrix result = weights.similarityT(A);
00810   result.invert(ierr);
00811   
00812   // blithely assuming the inverting never fails...
00813   return result;
00814 }
00815 
00816 void CSCSegAlgoSK::flipErrors( AlgebraicSymMatrix& a ) const {
00817 
00818   // The CSCSegment needs the error matrix re-arranged
00819 
00820   AlgebraicSymMatrix hold( a );
00821 
00822   // errors on slopes into upper left
00823   a(1,1) = hold(3,3);
00824   a(1,2) = hold(3,4);
00825   a(2,1) = hold(4,3);
00826   a(2,2) = hold(4,4);
00827   
00828   // errors on positions into lower right
00829   a(3,3) = hold(1,1);
00830   a(3,4) = hold(1,2);
00831   a(4,3) = hold(2,1);
00832   a(4,4) = hold(2,2);
00833 
00834   // off-diagonal elements remain unchanged
00835 
00836 }