CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoLocalMuon/CSCSegment/src/CSCSegAlgoTC.cc

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