CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/RecoLocalMuon/CSCSegment/src/CSCSegAlgoDF.cc

Go to the documentation of this file.
00001 
00007 #include "CSCSegAlgoDF.h"
00008 
00009 #include "Geometry/CSCGeometry/interface/CSCLayer.h"
00010 // For clhep Matrix::solve
00011 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
00012 
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 "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00019 
00020 #include "CSCSegAlgoPreClustering.h"
00021 #include "CSCSegAlgoShowering.h"
00022 
00023 #include <algorithm>
00024 #include <cmath>
00025 #include <iostream>
00026 #include <string>
00027 
00028 
00029 /* Constructor
00030  *
00031  */
00032 CSCSegAlgoDF::CSCSegAlgoDF(const edm::ParameterSet& ps) : CSCSegmentAlgorithm(ps), myName("CSCSegAlgoDF") {
00033         
00034   debug                  = ps.getUntrackedParameter<bool>("CSCSegmentDebug");
00035   minLayersApart         = ps.getParameter<int>("minLayersApart");
00036   minHitsPerSegment      = ps.getParameter<int>("minHitsPerSegment");
00037   dRPhiFineMax           = ps.getParameter<double>("dRPhiFineMax");
00038   dPhiFineMax            = ps.getParameter<double>("dPhiFineMax");
00039   tanThetaMax            = ps.getParameter<double>("tanThetaMax");
00040   tanPhiMax              = ps.getParameter<double>("tanPhiMax");        
00041   chi2Max                = ps.getParameter<double>("chi2Max");  
00042   preClustering          = ps.getUntrackedParameter<bool>("preClustering");
00043   minHitsForPreClustering= ps.getParameter<int>("minHitsForPreClustering");
00044   nHitsPerClusterIsShower= ps.getParameter<int>("nHitsPerClusterIsShower");
00045   Pruning                = ps.getUntrackedParameter<bool>("Pruning");
00046   maxRatioResidual       = ps.getParameter<double>("maxRatioResidualPrune");
00047 
00048   preCluster_            = new CSCSegAlgoPreClustering( ps );
00049   showering_             = new CSCSegAlgoShowering( ps );
00050 }
00051 
00052 
00053 /* Destructor
00054  *
00055  */
00056 CSCSegAlgoDF::~CSCSegAlgoDF() {
00057   delete preCluster_;
00058   delete showering_;
00059 }
00060 
00061 
00062 /* run
00063  *
00064  */
00065 std::vector<CSCSegment> CSCSegAlgoDF::run(const CSCChamber* aChamber, ChamberHitContainer rechits) {
00066 
00067   // Store chamber info in temp memory
00068   theChamber = aChamber; 
00069 
00070   int nHits = rechits.size();
00071 
00072   // Segments prior to pruning
00073   std::vector<CSCSegment> segments_temp;  
00074 
00075   if ( preClustering && nHits > minHitsForPreClustering ) {
00076     // This is where the segment origin is in the chamber on avg.
00077     std::vector<CSCSegment> testSegments;
00078     std::vector<ChamberHitContainer> clusteredHits = preCluster_->clusterHits(theChamber, rechits);
00079     // loop over the found clusters:
00080     for (std::vector<ChamberHitContainer>::iterator subrechits = clusteredHits.begin(); subrechits != clusteredHits.end(); ++subrechits ) {
00081       // build the subset of segments:
00082       std::vector<CSCSegment> segs = buildSegments( (*subrechits) );
00083       // add the found subset of segments to the collection of all segments in this chamber:
00084       segments_temp.insert( segments_temp.end(), segs.begin(), segs.end() );
00085     }
00086   } else {
00087     std::vector<CSCSegment> segs = buildSegments( rechits );
00088     // add the found subset of segments to the collection of all segments in this chamber:
00089     segments_temp.insert( segments_temp.end(), segs.begin(), segs.end() ); 
00090   }
00091 
00092   return segments_temp; 
00093 }
00094 
00095 
00096 /* This builds segments by first creating proto-segments from at least 3 hits.
00097  * We intend to try all possible pairs of hits to start segment building. 'All possible' 
00098  * means each hit lies on different layers in the chamber.  Once a hit has been assigned 
00099  * to a segment, we don't consider it again, THAT IS, FOR THE FIRST PASS ONLY !
00100  * In fact, this is one of the possible flaw with the SK algorithms as it sometimes manages
00101  * to build segments with the wrong starting points.  In the DF algorithm, the endpoints
00102  * are tested as the best starting points in a 2nd loop.
00103  *
00104  * Also, only a certain muonsPerChamberMax maximum number of segments can be produced in the chamber
00105  */
00106 std::vector<CSCSegment> CSCSegAlgoDF::buildSegments(ChamberHitContainer rechits) {
00107 
00108   // Clear buffer for segment vector
00109   std::vector<CSCSegment> segmentInChamber;
00110   segmentInChamber.clear();
00111 
00112   unsigned nHitInChamber = rechits.size();
00113   if ( nHitInChamber < 3 ) return segmentInChamber;
00114 
00115   LayerIndex layerIndex( nHitInChamber );
00116 
00117   unsigned nLayers = 0;
00118   int old_layer = -1;   
00119   for ( unsigned int i = 0; i < nHitInChamber; i++ ) {    
00120     int this_layer = rechits[i]->cscDetId().layer();
00121     layerIndex[i] = this_layer;
00122     if ( this_layer != old_layer ) {
00123       old_layer = this_layer;
00124       nLayers++;   
00125     }
00126   }
00127   
00128   if ( nLayers < 3 ) return segmentInChamber;
00129 
00130   double z1 = theChamber->layer(1)->position().z();
00131   double z6 = theChamber->layer(6)->position().z();
00132   
00133   if ( z1 > 0. ) {
00134     if ( z1 > z6 ) { 
00135       reverse( layerIndex.begin(), layerIndex.end() );
00136       reverse( rechits.begin(),    rechits.end() );
00137     }    
00138   }
00139   else if ( z1 < 0. ) {
00140     if ( z1 < z6 ) {
00141       reverse( layerIndex.begin(), layerIndex.end() );
00142       reverse( rechits.begin(),    rechits.end() );
00143     }    
00144   }
00145 
00146   // Showering muon
00147   if ( preClustering && int(nHitInChamber) > nHitsPerClusterIsShower && nLayers > 2 ) {
00148     CSCSegment segShower = showering_->showerSeg(theChamber, rechits);
00149 
00150     // Make sure have at least 3 hits...
00151     if ( segShower.nRecHits() < 3 ) return segmentInChamber;
00152 
00153     segmentInChamber.push_back(segShower);
00154 
00155     return segmentInChamber;
00156   }
00157 
00158 
00159   // Initialize flags that a given hit has been allocated to a segment
00160   BoolContainer used_ini(rechits.size(), false);
00161   usedHits = used_ini;
00162   
00163   ChamberHitContainerCIt ib = rechits.begin();
00164   ChamberHitContainerCIt ie = rechits.end();
00165     
00166   // Now Loop over hits within the chamber to find 1st seed for segment building
00167   for ( ChamberHitContainerCIt i1 = ib; i1 < ie; ++i1 ) {
00168     if ( usedHits[i1-ib] ) continue;
00169 
00170     const CSCRecHit2D* h1 = *i1;
00171     int layer1 = layerIndex[i1-ib];
00172     const CSCLayer* l1 = theChamber->layer(layer1);
00173     GlobalPoint gp1 = l1->toGlobal(h1->localPosition());
00174     LocalPoint lp1 = theChamber->toLocal(gp1);  
00175            
00176     // Loop over hits backward to find 2nd seed for segment building
00177     for ( ChamberHitContainerCIt i2 = ie-1; i2 > ib; --i2 ) {   
00178 
00179       if ( usedHits[i2-ib] ) continue;   // Hit has been used already
00180 
00181       int layer2 = layerIndex[i2-ib];
00182       if ( (layer2 - layer1) < minLayersApart ) continue;
00183 
00184       const CSCRecHit2D* h2 = *i2;
00185       const CSCLayer* l2 = theChamber->layer(layer2);
00186       GlobalPoint gp2 = l2->toGlobal(h2->localPosition());
00187       LocalPoint lp2 = theChamber->toLocal(gp2);  
00188 
00189       // Clear proto segment so it can be (re)-filled 
00190       protoSegment.clear();
00191 
00192       // localPosition is position of hit wrt layer (so local z = 0)
00193       protoIntercept = h1->localPosition();
00194 
00195       // We want hit wrt chamber (and local z will be != 0)
00196       float dz = gp2.z()-gp1.z();
00197       protoSlope_u = (lp2.x() - lp1.x())/dz ;
00198       protoSlope_v = (lp2.y() - lp1.y())/dz ;    
00199 
00200       // Test if entrance angle is roughly pointing towards IP
00201       if (fabs(protoSlope_v) > tanThetaMax) continue;
00202       if (fabs(protoSlope_u) > tanPhiMax ) continue;
00203      
00204       protoSegment.push_back(h1);
00205       protoSegment.push_back(h2);
00206         
00207       // Try adding hits to proto segment
00208       tryAddingHitsToSegment(rechits, i1, i2, layerIndex); 
00209         
00210       // Check no. of hits on segment to see if segment is large enough
00211       bool segok = true;
00212       unsigned iadd = 0;
00213 
00214       if (protoSegment.size() < minHitsPerSegment+iadd) segok = false;
00215   
00216       if ( Pruning && segok ) pruneFromResidual();
00217 
00218       // Check if segment satisfies chi2 requirement
00219       if (protoChi2 > chi2Max) segok = false;
00220 
00221       if ( segok ) {
00222 
00223         // Fill segment properties
00224        
00225         // Local direction
00226         double dz   = 1./sqrt(1. + protoSlope_u*protoSlope_u + protoSlope_v*protoSlope_v);
00227         double dx   = dz * protoSlope_u;
00228         double dy   = dz * protoSlope_v;
00229         LocalVector localDir(dx,dy,dz);
00230       
00231         // localDir may need sign flip to ensure it points outward from IP  
00232         double globalZpos    = ( theChamber->toGlobal( protoIntercept ) ).z();
00233         double globalZdir    = ( theChamber->toGlobal( localDir ) ).z();
00234         double directionSign = globalZpos * globalZdir;
00235         protoDirection       = (directionSign * localDir).unit();
00236 
00237         // Error matrix
00238         AlgebraicSymMatrix protoErrors = calculateError();     
00239         // but reorder components to match what's required by TrackingRecHit interface
00240         // i.e. slopes first, then positions
00241         flipErrors( protoErrors );
00242 
00243         CSCSegment temp(protoSegment, protoIntercept, protoDirection, protoErrors, protoChi2); 
00244 
00245         segmentInChamber.push_back(temp); 
00246 
00247         if (nHitInChamber-protoSegment.size() < 3) return segmentInChamber; 
00248         if (segmentInChamber.size() > 4) return segmentInChamber;
00249 
00250         // Flag used hits
00251         flagHitsAsUsed(rechits);
00252       } 
00253     } 
00254   }
00255   return segmentInChamber;
00256 }
00257 
00258 
00259 /* Method tryAddingHitsToSegment
00260  *
00261  * Look at left over hits and try to add them to proto segment by looking how far they
00262  * are from the segment in terms of the hit error matrix (so how many sigmas away).
00263  *
00264  */
00265 void CSCSegAlgoDF::tryAddingHitsToSegment( const ChamberHitContainer& rechits, 
00266                                            const ChamberHitContainerCIt i1, 
00267                                            const ChamberHitContainerCIt i2,
00268                                            LayerIndex layerIndex ) {
00269   
00270 /* Iterate over the layers with hits in the chamber
00271  * Skip the layers containing the segment endpoints on first pass, but then
00272  * try hits on layer containing the segment starting points on 2nd pass
00273  * if segment has >2 hits.  Once a hit is added to a layer, don't replace it 
00274  * until 2nd iteration
00275  */  
00276   
00277   ChamberHitContainerCIt ib = rechits.begin();
00278   ChamberHitContainerCIt ie = rechits.end();
00279   closeHits.clear();
00280 
00281   for ( ChamberHitContainerCIt i = ib; i != ie; ++i ) {
00282 
00283     if (i == i1 || i == i2 ) continue;
00284     if ( usedHits[i-ib] ) continue;   // Don't use hits already part of a segment.
00285 
00286     const CSCRecHit2D* h = *i;
00287     int layer = layerIndex[i-ib];
00288     int layer1 = layerIndex[i1-ib];
00289     int layer2 = layerIndex[i2-ib];
00290 
00291     // Low multiplicity case
00292     if (rechits.size() < 9) {
00293       if ( isHitNearSegment( h ) ) {
00294         if ( !hasHitOnLayer(layer) ) {
00295           addHit(h, layer);
00296         } else {
00297           closeHits.push_back(h);
00298         }
00299       }
00300 
00301     // High multiplicity case
00302     } else { 
00303       if ( isHitNearSegment( h ) ) {
00304         if ( !hasHitOnLayer(layer) ) {
00305           addHit(h, layer);
00306           updateParameters();
00307         // Don't change the starting points at this stage !!!
00308         } else {
00309           closeHits.push_back(h);
00310           if (layer != layer1 && layer != layer2 ) compareProtoSegment(h, layer);
00311         }
00312       }
00313     }
00314   }
00315  
00316   if ( int(protoSegment.size()) < 3) return;
00317 
00318   updateParameters();
00319 
00320   // 2nd pass to remove biases 
00321   // This time, also consider changing the endpoints
00322   for ( ChamberHitContainerCIt i = closeHits.begin() ; i != closeHits.end(); ++i ) {      
00323     const CSCRecHit2D* h = *i;      
00324     int layer = (*i)->cscDetId().layer();     
00325     compareProtoSegment(h, layer); 
00326   } 
00327 
00328 }
00329 
00330 
00331 /* isHitNearSegment
00332  *
00333  * Compare rechit with expected position from proto_segment
00334  */
00335 bool CSCSegAlgoDF::isHitNearSegment( const CSCRecHit2D* hit) const {
00336 
00337   const CSCLayer* layer = theChamber->layer(hit->cscDetId().layer());
00338 
00339   // hit phi position in global coordinates
00340   GlobalPoint Hgp = layer->toGlobal(hit->localPosition());
00341   double Hphi = Hgp.phi();                                
00342   if (Hphi < 0.) Hphi += 2.*M_PI;
00343   LocalPoint Hlp = theChamber->toLocal(Hgp);
00344   double z = Hlp.z();  
00345 
00346   double LocalX = protoIntercept.x() + protoSlope_u * z;
00347   double LocalY = protoIntercept.y() + protoSlope_v * z;
00348   LocalPoint Slp(LocalX, LocalY, z);
00349   GlobalPoint Sgp = theChamber->toGlobal(Slp); 
00350   double Sphi = Sgp.phi();
00351   if (Sphi < 0.) Sphi += 2.*M_PI;
00352   double R = sqrt(Sgp.x()*Sgp.x() + Sgp.y()*Sgp.y());
00353   
00354   double deltaPhi = Sphi - Hphi;
00355   if (deltaPhi >  2.*M_PI) deltaPhi -= 2.*M_PI;
00356   if (deltaPhi < -2.*M_PI) deltaPhi += 2.*M_PI;
00357   if (deltaPhi < 0.) deltaPhi = -deltaPhi; 
00358 
00359   double RdeltaPhi = R * deltaPhi;
00360 
00361   if (RdeltaPhi < dRPhiFineMax && deltaPhi < dPhiFineMax ) return true;
00362 
00363   return false;
00364 }
00365 
00366 
00367 /* Method addHit
00368  *
00369  * Test if can add hit to proto segment. If so, try to add it.
00370  *
00371  */
00372 bool CSCSegAlgoDF::addHit(const CSCRecHit2D* aHit, int layer) {
00373   
00374   // Return true if hit was added successfully and then parameters are updated.
00375   // Return false if there is already a hit on the same layer, or insert failed.
00376   
00377   bool ok = true;
00378   
00379   // Test that we are not trying to add the same hit again
00380   for ( ChamberHitContainer::const_iterator it = protoSegment.begin(); it != protoSegment.end(); it++ ) 
00381     if ( aHit == (*it)  ) return false;
00382   
00383   protoSegment.push_back(aHit);
00384 
00385   return ok;
00386 }    
00387 
00388 
00389 /* Method updateParameters
00390  *      
00391  * Perform a simple Least Square Fit on proto segment to determine slope and intercept
00392  *
00393  */   
00394 void CSCSegAlgoDF::updateParameters() {
00395 
00396   // Compute slope from Least Square Fit    
00397   CLHEP::HepMatrix M(4,4,0);
00398   CLHEP::HepVector B(4,0);
00399 
00400   ChamberHitContainer::const_iterator ih;
00401   
00402   for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
00403     
00404     const CSCRecHit2D& hit = (**ih);
00405     const CSCLayer* layer  = theChamber->layer(hit.cscDetId().layer());
00406     GlobalPoint gp         = layer->toGlobal(hit.localPosition());
00407     LocalPoint  lp         = theChamber->toLocal(gp); 
00408     
00409     double u = lp.x();
00410     double v = lp.y();
00411     double z = lp.z();
00412     
00413     // ptc: Covariance matrix of local errors 
00414     CLHEP::HepMatrix IC(2,2);
00415     IC(1,1) = hit.localPositionError().xx();
00416     IC(1,2) = hit.localPositionError().xy();
00417     IC(2,2) = hit.localPositionError().yy();
00418     IC(2,1) = IC(1,2); // since Cov is symmetric
00419     
00420     // ptc: Invert covariance matrix (and trap if it fails!)
00421     int ierr = 0;
00422     IC.invert(ierr); // inverts in place
00423     if (ierr != 0) {
00424       LogDebug("CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC << "\n";      
00425     }
00426     
00427     M(1,1) += IC(1,1);
00428     M(1,2) += IC(1,2);
00429     M(1,3) += IC(1,1) * z;
00430     M(1,4) += IC(1,2) * z;
00431     B(1)   += u * IC(1,1) + v * IC(1,2);
00432     
00433     M(2,1) += IC(2,1);
00434     M(2,2) += IC(2,2);
00435     M(2,3) += IC(2,1) * z;
00436     M(2,4) += IC(2,2) * z;
00437     B(2)   += u * IC(2,1) + v * IC(2,2);
00438     
00439     M(3,1) += IC(1,1) * z;
00440     M(3,2) += IC(1,2) * z;
00441     M(3,3) += IC(1,1) * z * z;
00442     M(3,4) += IC(1,2) * z * z;
00443     B(3)   += ( u * IC(1,1) + v * IC(1,2) ) * z;
00444     
00445     M(4,1) += IC(2,1) * z;
00446     M(4,2) += IC(2,2) * z;
00447     M(4,3) += IC(2,1) * z * z;
00448     M(4,4) += IC(2,2) * z * z;
00449     B(4)   += ( u * IC(2,1) + v * IC(2,2) ) * z;
00450   }
00451   
00452   CLHEP::HepVector p = solve(M, B);
00453   
00454   // Update member variables 
00455   // Note that origin has local z = 0
00456 
00457   protoIntercept = LocalPoint(p(1), p(2), 0.);
00458   protoSlope_u = p(3);
00459   protoSlope_v = p(4);
00460 
00461 
00462   // Determine Chi^2 for the proto wire segment
00463   
00464   double chsq = 0.;
00465   
00466   for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
00467     
00468     const CSCRecHit2D& hit = (**ih);
00469     const CSCLayer* layer  = theChamber->layer(hit.cscDetId().layer());
00470     GlobalPoint gp         = layer->toGlobal(hit.localPosition());
00471     LocalPoint lp          = theChamber->toLocal(gp);
00472     
00473     double u = lp.x();
00474     double v = lp.y();
00475     double z = lp.z();
00476     
00477     double du = protoIntercept.x() + protoSlope_u * z - u;
00478     double dv = protoIntercept.y() + protoSlope_v * z - v;
00479     
00480     CLHEP::HepMatrix IC(2,2);
00481     IC(1,1) = hit.localPositionError().xx();
00482     IC(1,2) = hit.localPositionError().xy();
00483     IC(2,2) = hit.localPositionError().yy();
00484     IC(2,1) = IC(1,2);
00485     
00486     // Invert covariance matrix
00487     int ierr = 0;
00488     IC.invert(ierr);
00489     if (ierr != 0) {
00490       LogDebug("CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC << "\n";      
00491     }
00492     chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
00493   }
00494   protoChi2 = chsq;
00495 }
00496 
00497 
00498 /* hasHitOnLayer
00499  *
00500  * Just make sure hit to be added to layer comes from different layer than those in proto segment   
00501  */
00502 bool CSCSegAlgoDF::hasHitOnLayer(int layer) const {
00503   
00504   // Is there already a hit on this layer?
00505   for ( ChamberHitContainerCIt it = protoSegment.begin(); it != protoSegment.end(); it++ )
00506     if ( (*it)->cscDetId().layer() == layer ) return true;
00507   
00508   return false;
00509 }
00510 
00511 
00512 /* Method compareProtoSegment
00513  *      
00514  * For hit coming from the same layer of an existing hit within the proto segment
00515  * test if achieve better chi^2 by using this hit than the other
00516  *
00517  */ 
00518 void CSCSegAlgoDF::compareProtoSegment(const CSCRecHit2D* h, int layer) {
00519   
00520   // Store old segment first
00521   double old_protoChi2                  = protoChi2;
00522   LocalPoint old_protoIntercept         = protoIntercept;
00523   float old_protoSlope_u                = protoSlope_u;
00524   float old_protoSlope_v                = protoSlope_v;
00525   LocalVector old_protoDirection        = protoDirection;
00526   ChamberHitContainer old_protoSegment  = protoSegment;
00527  
00528 
00529   // Try adding the hit to existing segment, and remove old one existing in same layer
00530   ChamberHitContainer::iterator it;
00531   for ( it = protoSegment.begin(); it != protoSegment.end(); ) {
00532     if ( (*it)->cscDetId().layer() == layer ) {
00533       it = protoSegment.erase(it);
00534     } else {
00535       ++it;
00536     }
00537   }
00538   bool ok = addHit(h, layer);
00539 
00540   if (ok) updateParameters();
00541   
00542   if ( (protoChi2 > old_protoChi2) || ( !ok ) ) {
00543     protoChi2       = old_protoChi2;
00544     protoIntercept  = old_protoIntercept;
00545     protoSlope_u    = old_protoSlope_u;
00546     protoSlope_v    = old_protoSlope_v;
00547     protoDirection  = old_protoDirection;
00548     protoSegment    = old_protoSegment;
00549   }
00550 }
00551 
00552 
00553 /* Method flagHitsAsUsed
00554  *
00555  * Flag hits which have entered segment building so we don't reuse them.
00556  * Also flag does which were very close to segment to reduce combinatorics
00557  */
00558 void CSCSegAlgoDF::flagHitsAsUsed(const ChamberHitContainer& rechitsInChamber) {
00559   
00560   // Flag hits on segment as used
00561   ChamberHitContainerCIt ib = rechitsInChamber.begin();
00562   ChamberHitContainerCIt hi, iu;
00563   
00564   for ( hi = protoSegment.begin(); hi != protoSegment.end(); ++hi ) {
00565     for ( iu = ib; iu != rechitsInChamber.end(); ++iu ) {
00566       if (*hi == *iu) usedHits[iu-ib] = true;
00567     }
00568   }
00569   // Don't reject hits marked as "nearby" for now.
00570   // So this is bypassed at all times for now !!!
00571   // Perhaps add back to speed up algorithm some more
00572   if (closeHits.size() > 0) return;  
00573   for ( hi = closeHits.begin(); hi != closeHits.end(); ++hi ) {
00574     for ( iu = ib; iu != rechitsInChamber.end(); ++iu ) {
00575       if (*hi == *iu) usedHits[iu-ib] = true;
00576     }
00577   }
00578 
00579 }
00580 
00581 
00582 AlgebraicSymMatrix CSCSegAlgoDF::weightMatrix() const {
00583 
00584   std::vector<const CSCRecHit2D*>::const_iterator it;
00585   int nhits = protoSegment.size();
00586   AlgebraicSymMatrix matrix(2*nhits, 0);
00587   int row = 0;
00588 
00589   for (it = protoSegment.begin(); it != protoSegment.end(); ++it) {
00590 
00591     const CSCRecHit2D& hit = (**it);
00592     ++row;
00593     matrix(row, row)   = hit.localPositionError().xx();
00594     matrix(row, row+1) = hit.localPositionError().xy();
00595     ++row;
00596     matrix(row, row-1) = hit.localPositionError().xy();
00597     matrix(row, row)   = hit.localPositionError().yy();
00598   }
00599   int ierr;
00600   matrix.invert(ierr);
00601   return matrix;
00602 }
00603 
00604 
00605 CLHEP::HepMatrix CSCSegAlgoDF::derivativeMatrix() const {
00606 
00607   ChamberHitContainer::const_iterator it;
00608   int nhits = protoSegment.size();
00609   CLHEP::HepMatrix matrix(2*nhits, 4);
00610   int row = 0;
00611 
00612   for(it = protoSegment.begin(); it != protoSegment.end(); ++it) {
00613 
00614     const CSCRecHit2D& hit = (**it);
00615     const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
00616     GlobalPoint gp = layer->toGlobal(hit.localPosition());
00617     LocalPoint lp = theChamber->toLocal(gp);
00618     float z = lp.z();
00619     ++row;
00620     matrix(row, 1) = 1.;
00621     matrix(row, 3) = z;
00622     ++row;
00623     matrix(row, 2) = 1.;
00624     matrix(row, 4) = z;
00625   }
00626   return matrix;
00627 }
00628 
00629 /* calculateError
00630  *
00631  */
00632 AlgebraicSymMatrix CSCSegAlgoDF::calculateError() const {
00633 
00634   AlgebraicSymMatrix weights = weightMatrix();
00635   AlgebraicMatrix A = derivativeMatrix();
00636 
00637   // (AT W A)^-1
00638   // from http://www.phys.ufl.edu/~avery/fitting.html, part I
00639   int ierr;
00640   AlgebraicSymMatrix result = weights.similarityT(A);
00641   result.invert(ierr);
00642 
00643   // blithely assuming the inverting never fails...
00644   return result;
00645 }
00646 
00647 void CSCSegAlgoDF::flipErrors( AlgebraicSymMatrix& a ) const {
00648 
00649   // The CSCSegment needs the error matrix re-arranged
00650 
00651   AlgebraicSymMatrix hold( a );
00652 
00653   // errors on slopes into upper left
00654   a(1,1) = hold(3,3);
00655   a(1,2) = hold(3,4);
00656   a(2,1) = hold(4,3);
00657   a(2,2) = hold(4,4);
00658 
00659   // errors on positions into lower right
00660   a(3,3) = hold(1,1);
00661   a(3,4) = hold(1,2);
00662   a(4,3) = hold(2,1);
00663   a(4,4) = hold(2,2);
00664 
00665   // off-diagonal elements remain unchanged
00666 
00667 }
00668 
00669 
00670 // Try to clean up segments by quickly looking at residuals
00671 void CSCSegAlgoDF::pruneFromResidual(){
00672 
00673   // Only prune if have at least 5 hits 
00674   if ( protoSegment.size() < 5 ) return ;
00675 
00676 
00677   // Now Study residuals
00678       
00679   float maxResidual = 0.;
00680   float sumResidual = 0.;
00681   int nHits = 0;
00682   int badIndex = -1;
00683   int j = 0;
00684 
00685 
00686   ChamberHitContainer::const_iterator ih;
00687 
00688   for ( ih = protoSegment.begin(); ih != protoSegment.end(); ++ih ) {
00689     const CSCRecHit2D& hit = (**ih);
00690     const CSCLayer* layer  = theChamber->layer(hit.cscDetId().layer());
00691     GlobalPoint gp         = layer->toGlobal(hit.localPosition());
00692     LocalPoint lp          = theChamber->toLocal(gp);
00693 
00694     double u = lp.x();
00695     double v = lp.y();
00696     double z = lp.z();
00697 
00698     double du = protoIntercept.x() + protoSlope_u * z - u;
00699     double dv = protoIntercept.y() + protoSlope_v * z - v;
00700 
00701     float residual = sqrt(du*du + dv*dv);
00702 
00703     sumResidual += residual;
00704     nHits++;
00705     if ( residual > maxResidual ) {
00706       maxResidual = residual;
00707       badIndex = j;
00708     }
00709     j++;
00710   }
00711 
00712   float corrAvgResidual = (sumResidual - maxResidual)/(nHits -1);
00713 
00714   // Keep all hits 
00715   if ( maxResidual/corrAvgResidual < maxRatioResidual ) return;
00716 
00717 
00718   // Drop worse hit and recompute segment properties + fill
00719 
00720   ChamberHitContainer newProtoSegment;
00721 
00722   j = 0;
00723   for ( ih = protoSegment.begin(); ih != protoSegment.end(); ++ih ) {
00724     if ( j != badIndex ) newProtoSegment.push_back(*ih);
00725     j++;
00726   }
00727   
00728   protoSegment.clear();
00729 
00730   for ( ih = newProtoSegment.begin(); ih != newProtoSegment.end(); ++ih ) {
00731     protoSegment.push_back(*ih);
00732   }
00733 
00734   // Update segment parameters
00735   updateParameters();
00736 
00737 }
00738 
00739 
00740 /*
00741  * Order the hits such that 2nd one is closest in x,y to first seed hit in global coordinates
00742  */
00743 void CSCSegAlgoDF::orderSecondSeed( GlobalPoint gp1,
00744                                            const ChamberHitContainerCIt i1, 
00745                                            const ChamberHitContainerCIt i2, 
00746                                            const ChamberHitContainer& rechits, 
00747                                            LayerIndex layerIndex ) {
00748 
00749   secondSeedHits.clear();
00750 
00751   ChamberHitContainerCIt ib = rechits.begin();
00752   ChamberHitContainerCIt ie = rechits.end();
00753 
00754   //  int layer1 = layerIndex[i1-ib];
00755   //  int layer2 = layerIndex[i2-ib];
00756 
00757 
00758   // Now fill vector of rechits closest to center of mass:
00759   // secondSeedHitsIdx.clear() = 0;
00760 
00761   // Loop over all hits and find hit closest to 1st seed.
00762   for ( ChamberHitContainerCIt i2 = ie-1; i2 > i1; --i2 ) {     
00763 
00764 
00765   }
00766 
00767         
00768 }