CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoEgamma/EgammaElectronAlgos/src/SiStripElectronAlgo.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     EgammaElectronAlgos
00004 // Class  :     SiStripElectronAlgo
00005 // 
00006 // Implementation:
00007 //     <Notes on implementation>
00008 //
00009 // Original Author:  Jim Pivarski
00010 //         Created:  Fri May 26 16:12:04 EDT 2006
00011 // $Id: SiStripElectronAlgo.cc,v 1.31 2010/09/21 17:06:14 chamont Exp $
00012 //
00013 
00014 // system include files
00015 
00016 // user include files
00017 #include "DataFormats/EgammaCandidates/interface/SiStripElectronFwd.h"
00018 #include "RecoEgamma/EgammaElectronAlgos/interface/SiStripElectronAlgo.h"
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 
00021 #include "Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h"
00022 #include "Geometry/CommonTopologies/interface/RectangularStripTopology.h"
00023 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00024 
00025 #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
00026 #include "RecoTracker/TrackProducer/interface/TrackingRecHitLessFromGlobalPosition.h"
00027 
00028 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00029 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00030 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00031 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00032 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00033 
00034 
00035 //
00036 // constants, enums and typedefs
00037 //
00038 
00039 //
00040 // static data member definitions
00041 //
00042 namespace {
00043   struct CompareHits {
00044     CompareHits( const TrackingRecHitLessFromGlobalPosition& iLess ) :
00045       less_(iLess) {}
00046     bool operator()(const TrackingRecHit* iLHS, const TrackingRecHit* iRHS) {
00047       return less_(*iLHS,*iRHS);
00048     }
00049          
00050     TrackingRecHitLessFromGlobalPosition less_;
00051   };
00052 }
00053 
00054 
00055 //
00056 // constructors and destructor
00057 //
00058 SiStripElectronAlgo::SiStripElectronAlgo(unsigned int maxHitsOnDetId,
00059                                          double originUncertainty,
00060                                          double phiBandWidth,
00061                                          double maxNormResid,
00062                                          unsigned int minHits,
00063                                          double maxReducedChi2)
00064   : maxHitsOnDetId_(maxHitsOnDetId)
00065   , originUncertainty_(originUncertainty)
00066   , phiBandWidth_(phiBandWidth)
00067   , maxNormResid_(maxNormResid)
00068   , minHits_(minHits)
00069   , maxReducedChi2_(maxReducedChi2)
00070 {
00071 }
00072 
00073 // SiStripElectronAlgo::SiStripElectronAlgo(const SiStripElectronAlgo& rhs)
00074 // {
00075 //    // do actual copying here;
00076 // }
00077 
00078 SiStripElectronAlgo::~SiStripElectronAlgo()
00079 {
00080 }
00081 
00082 //
00083 // assignment operators
00084 //
00085 // const SiStripElectronAlgo& SiStripElectronAlgo::operator=(const SiStripElectronAlgo& rhs)
00086 // {
00087 //   //An exception safe implementation is
00088 //   SiStripElectronAlgo temp(rhs);
00089 //   swap(rhs);
00090 //
00091 //   return *this;
00092 // }
00093 
00094 //
00095 // member functions
00096 //
00097 
00098 void SiStripElectronAlgo::prepareEvent(const edm::ESHandle<TrackerGeometry>& tracker,
00099                                        const edm::Handle<SiStripRecHit2DCollection>& rphiHits,
00100                                        const edm::Handle<SiStripRecHit2DCollection>& stereoHits,
00101                                        const edm::Handle<SiStripMatchedRecHit2DCollection>& matchedHits,
00102                                        const edm::ESHandle<MagneticField>& magneticField)
00103 {
00104   LogDebug("") << " In prepareEvent " ; 
00105 
00106   tracker_p_ = tracker.product();
00107   rphiHits_p_ = rphiHits.product();
00108   stereoHits_p_ = stereoHits.product();
00109   matchedHits_p_ = matchedHits.product();
00110   magneticField_p_ = magneticField.product();
00111 
00112   rphiHits_hp_ = &rphiHits;
00113   stereoHits_hp_ = &stereoHits;
00114 
00115   // Keep a table that relates hit pointers to their index (key) in the collections
00116   rphiKey_.clear();
00117   stereoKey_.clear();
00118   // Keep track of which hits have been used already (so a hit is assigned to only one electron)
00119   hitUsed_.clear();
00120   matchedHitUsed_.clear();
00121 
00122   unsigned int counter = 0;
00123   for (SiStripRecHit2DCollection::DataContainer::const_iterator it = rphiHits_p_->data().begin();  it != rphiHits_p_->data().end();  ++it) {
00124     rphiKey_[&(*it)] = counter;
00125     hitUsed_[&(*it)] = false;
00126     counter++;
00127   }
00128 
00129   counter = 0;
00130   for (SiStripRecHit2DCollection::DataContainer::const_iterator it = stereoHits_p_->data().begin();  it != stereoHits_p_->data().end();  ++it) {
00131     stereoKey_[&(*it)] = counter;
00132     hitUsed_[&(*it)] = false;
00133     counter++;
00134   }
00135 
00136   counter = 0;
00137   for (SiStripMatchedRecHit2DCollection::DataContainer::const_iterator it = matchedHits_p_->data().begin();  it != matchedHits_p_->data().end();  ++it) {
00138     matchedKey_[&(*it)] = counter;
00139     matchedHitUsed_[&(*it)] = false;
00140     counter++;
00141   }
00142   
00143   LogDebug("") << " Leaving prepareEvent " ;
00144 }
00145 
00146 
00147 
00148 // returns true iff an electron was found
00149 // inserts electrons and trackcandiates into electronOut and trackCandidateOut
00150 bool SiStripElectronAlgo::findElectron(reco::SiStripElectronCollection& electronOut,
00151                                        TrackCandidateCollection& trackCandidateOut,
00152                                        const reco::SuperClusterRef& superclusterIn)
00153 {
00154   // Try each of the two charge hypotheses, but only take one
00155   bool electronSuccess = projectPhiBand(-1., superclusterIn);
00156   bool positronSuccess = projectPhiBand( 1., superclusterIn);
00157 
00158   // electron hypothesis did better than electron
00159   if ((electronSuccess  &&  !positronSuccess)  ||
00160       (electronSuccess  &&  positronSuccess  &&  redchi2_neg_ <= redchi2_pos_)) {
00161       
00162     // Initial uncertainty for tracking
00163     AlgebraicSymMatrix errors(5,1);  // makes identity 5x5 matrix, indexed from (1,1) to (5,5)
00164     errors(1,1) = 3.;      // uncertainty**2 in 1/momentum
00165     errors(2,2) = 0.01;    // uncertainty**2 in lambda (lambda == pi/2 - polar angle theta)
00166     errors(3,3) = 0.0001;  // uncertainty**2 in phi
00167     errors(4,4) = 0.01;    // uncertainty**2 in x_transverse (where x is in cm)
00168     errors(5,5) = 0.01;    // uncertainty**2 in y_transverse (where y is in cm)
00169 
00170     // JED Debugging possible double hit sort problem
00171     std::ostringstream debugstr6;
00172     debugstr6 << " HERE BEFORE SORT electron case " << " \n" ;
00173     for (std::vector<const TrackingRecHit*>::iterator itHit=outputHits_neg_.begin(); 
00174          itHit != outputHits_neg_.end(); ++itHit) {
00175       debugstr6 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
00176                 <<"    Local Position: "<<(*itHit)->localPosition()<<" \n"
00177                 <<"    Global Rho:  "
00178                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
00179                 <<" Phi "
00180                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
00181                 << " Z "
00182                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";
00183     }
00184     // end of dump 
00185     
00186     
00187     
00188     // JED call dump 
00189     debugstr6 << " Calling sort alongMomentum " << " \n"; 
00190     
00191     std::sort(outputHits_neg_.begin(), outputHits_neg_.end(),
00192               CompareHits(TrackingRecHitLessFromGlobalPosition(((TrackingGeometry*)(tracker_p_)), alongMomentum)));
00193     
00194 
00195     debugstr6 << " Done with sort " << " \n";
00196     
00197     debugstr6 << " HERE AFTER SORT electron case " << " \n";
00198     for (std::vector<const TrackingRecHit*>::iterator itHit=outputHits_neg_.begin(); 
00199          itHit != outputHits_neg_.end(); ++itHit) {
00200       debugstr6 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
00201                 <<"    Local Position: "<<(*itHit)->localPosition()<<" \n"
00202                 <<"    Global Rho:  "
00203                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
00204                 <<" Phi "
00205                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
00206                 << " Z "
00207                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";; 
00208     }
00209     // end of dump 
00210 
00211     LogDebug("")<< debugstr6.str();
00212 
00213     
00214     //create an OwnVector needed by the classes which will be stored to the Event
00215     edm::OwnVector<TrackingRecHit> hits;
00216     for(std::vector<const TrackingRecHit*>::iterator itHit =outputHits_neg_.begin();
00217         itHit != outputHits_neg_.end();
00218         ++itHit) {
00219       hits.push_back( (*itHit)->clone());
00220       if( !(hitUsed_.find(*itHit) != hitUsed_.end()) ) {
00221         LogDebug("") << " Assert failure " ;
00222         assert(hitUsed_.find(*itHit) != hitUsed_.end());
00223       }
00224       hitUsed_[*itHit] = true;
00225     }
00226     
00227     TrajectoryStateOnSurface state(
00228                                    GlobalTrajectoryParameters(position_neg_, momentum_neg_, -1, magneticField_p_),
00229                                    CurvilinearTrajectoryError(errors),
00230                                    tracker_p_->idToDet(innerhit_neg_->geographicalId())->surface());
00231     
00232     TrajectoryStateTransform transformer;
00233     PTrajectoryStateOnDet* PTraj = transformer.persistentState(state, innerhit_neg_->geographicalId().rawId());
00234     TrajectorySeed trajectorySeed(*PTraj, hits, alongMomentum);
00235     trackCandidateOut.push_back(TrackCandidate(hits, trajectorySeed, *PTraj));
00236     
00237     electronOut.push_back(reco::SiStripElectron(superclusterIn,
00238                                                 -1,
00239                                                 outputRphiHits_neg_,
00240                                                 outputStereoHits_neg_,
00241                                                 phiVsRSlope_neg_,
00242                                                 slope_neg_,
00243                                                 intercept_neg_ + superclusterIn->position().phi(),
00244                                                 chi2_neg_,
00245                                                 ndof_neg_,
00246                                                 correct_pT_neg_,
00247                                                 pZ_neg_,
00248                                                 zVsRSlope_neg_,
00249                                                 numberOfStereoHits_neg_,
00250                                                 numberOfBarrelRphiHits_neg_,
00251                                                 numberOfEndcapZphiHits_neg_));
00252     
00253     delete PTraj; // ShR: fix memory leak reported per perfrmance task force
00254     
00255     return true;
00256   }
00257   
00258   // positron hypothesis did better than electron
00259   if ((!electronSuccess  &&  positronSuccess)  ||
00260       (electronSuccess  &&  positronSuccess  &&  redchi2_neg_ > redchi2_pos_)) {
00261       
00262     // Initial uncertainty for tracking
00263     AlgebraicSymMatrix errors(5,1);  // makes identity 5x5 matrix, indexed from (1,1) to (5,5)
00264     errors(1,1) = 3.;      // uncertainty**2 in 1/momentum
00265     errors(2,2) = 0.01;    // uncertainty**2 in lambda (lambda == pi/2 - polar angle theta)
00266     errors(3,3) = 0.0001;  // uncertainty**2 in phi
00267     errors(4,4) = 0.01;    // uncertainty**2 in x_transverse (where x is in cm)
00268     errors(5,5) = 0.01;    // uncertainty**2 in y_transverse (where y is in cm)
00269 
00270     // JED Debugging possible double hit sort problem
00271     std::ostringstream debugstr7;
00272     debugstr7 << " HERE BEFORE SORT Positron case " << " \n";
00273     for (std::vector<const TrackingRecHit*>::iterator itHit = outputHits_pos_.begin(); 
00274          itHit != outputHits_pos_.end(); ++itHit) {
00275       debugstr7 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
00276                 <<"    Local Position: "<<(*itHit)->localPosition()<<" \n"
00277                 <<"    Global Rho:  "
00278                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
00279                 <<" Phi "
00280                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
00281                 << " Z "
00282                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";
00283     }
00284     // end of dump 
00285     
00286     debugstr7 << " Calling sort alongMomentum " << " \n"; 
00287     
00288     std::sort(outputHits_pos_.begin(), outputHits_pos_.end(),
00289               CompareHits(TrackingRecHitLessFromGlobalPosition(((TrackingGeometry*)(tracker_p_)), alongMomentum)));
00290     
00291     debugstr7 << " Done with sort " << " \n";
00292     
00293     // JED Debugging possible double hit sort problem
00294     debugstr7 << " HERE AFTER SORT Positron case " << " \n";
00295     for (std::vector<const TrackingRecHit*>::iterator itHit = outputHits_pos_.begin(); 
00296          itHit != outputHits_pos_.end(); ++itHit) {
00297       debugstr7 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
00298                 <<"    Local Position: "<<(*itHit)->localPosition()<<" \n"
00299                 <<"    Global Rho:  "
00300                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
00301                 <<" Phi "
00302                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
00303                 << " Z "
00304                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n"; 
00305     }
00306     // end of dump 
00307     LogDebug("") << debugstr7.str();
00308 
00309     //create an OwnVector needed by the classes which will be stored to the Event
00310     edm::OwnVector<TrackingRecHit> hits;
00311     for(std::vector<const TrackingRecHit*>::iterator itHit =outputHits_pos_.begin();
00312         itHit != outputHits_pos_.end();
00313         ++itHit) {
00314       hits.push_back( (*itHit)->clone());
00315       assert(hitUsed_.find(*itHit) != hitUsed_.end());
00316       hitUsed_[*itHit] = true;
00317     }
00318     
00319     TrajectoryStateOnSurface state(
00320                                    GlobalTrajectoryParameters(position_pos_, momentum_pos_, 1, magneticField_p_),
00321                                    CurvilinearTrajectoryError(errors),
00322                                    tracker_p_->idToDet(innerhit_pos_->geographicalId())->surface());
00323     
00324     TrajectoryStateTransform transformer;
00325     PTrajectoryStateOnDet* PTraj = transformer.persistentState(state, innerhit_pos_->geographicalId().rawId());
00326     TrajectorySeed trajectorySeed(*PTraj, hits, alongMomentum);
00327     trackCandidateOut.push_back(TrackCandidate(hits, trajectorySeed, *PTraj));
00328     
00329     electronOut.push_back(reco::SiStripElectron(superclusterIn,
00330                                                 1,
00331                                                 outputRphiHits_pos_,
00332                                                 outputStereoHits_pos_,
00333                                                 phiVsRSlope_pos_,
00334                                                 slope_pos_,
00335                                                 intercept_pos_ + superclusterIn->position().phi(),
00336                                                 chi2_pos_,
00337                                                 ndof_pos_,
00338                                                 correct_pT_pos_,
00339                                                 pZ_pos_,
00340                                                 zVsRSlope_pos_,
00341                                                 numberOfStereoHits_pos_,
00342                                                 numberOfBarrelRphiHits_pos_,
00343                                                 numberOfEndcapZphiHits_pos_));
00344     
00345     delete PTraj; // JED: fix of 2nd memory leak reported per perfrmance task force
00346     
00347     return true;
00348   }
00349   
00350   return false;
00351 }
00352 
00353 // inserts pointers to good hits into hitPointersOut
00354 // selects hits on DetIds that have no more than maxHitsOnDetId_
00355 // selects from stereo if stereo == true, rphi otherwise
00356 // selects from TID or TEC if endcap == true, TIB or TOB otherwise
00357 void SiStripElectronAlgo::coarseHitSelection(std::vector<const SiStripRecHit2D*>& hitPointersOut,
00358                                              bool stereo, bool endcap)
00359 {
00360   // This function is not time-efficient.  If you want to improve the
00361   // speed of the algorithm, you'll probably want to change this
00362   // function.  There may be a more efficienct way to extract hits,
00363   // and it would definitely help to put a geographical cut on the
00364   // DetIds.  (How does one determine the global position of a given
00365   // DetId?  Is tracker_p_->idToDet(id)->surface().toGlobal(LocalPosition(0,0,0))
00366   // expensive?)
00367 
00368   // Loop over the detector ids
00369   SiStripRecHit2DCollection::const_iterator itdet = (stereo ? stereoHits_p_->begin() : rphiHits_p_->begin());
00370   SiStripRecHit2DCollection::const_iterator eddet = (stereo ? stereoHits_p_->end()   : rphiHits_p_->end()  );
00371   for (; itdet != eddet; ++itdet) {
00372     // Get the hits on this detector id
00373     SiStripRecHit2DCollection::DetSet hits = *itdet;
00374     DetId id(hits.detId());
00375 
00376     // Count the number of hits on this detector id
00377     unsigned int numberOfHits = hits.size();
00378       
00379     // Only take the hits if there aren't too many
00380     // (Would it be better to loop only once, fill a temporary list,
00381     // and copy that if numberOfHits <= maxHitsOnDetId_?)
00382     if (numberOfHits <= maxHitsOnDetId_) {
00383       for (SiStripRecHit2DCollection::DetSet::const_iterator hit = hits.begin();  
00384            hit != hits.end();  ++hit) {
00385         // check that hit is valid first !
00386         if(!(*hit).isValid()) {
00387           LogDebug("") << " InValid hit skipped in coarseHitSelection " << std::endl ;
00388           continue ;
00389         }
00390         std::string theDet = "null";
00391         int theLayer = -999;
00392         bool isStereoDet = false ;
00393         if(tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TIB) { 
00394           theDet = "TIB" ;
00395           theLayer = TIBDetId(id).layer(); 
00396           if(TIBDetId(id).stereo()==1) { isStereoDet = true ; }
00397         } else if
00398           (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TOB) { 
00399           theDet = "TOB" ;
00400           theLayer = TOBDetId(id).layer(); 
00401           if(TOBDetId(id).stereo()==1) { isStereoDet = true ; }
00402         }else if
00403           (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TID) { 
00404           theDet = "TID" ;
00405           theLayer = TIDDetId(id).wheel();  // or ring  ?
00406           if(TIDDetId(id).stereo()==1) { isStereoDet = true ; }
00407         }else if
00408           (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TEC) { 
00409           theDet = "TEC" ;
00410           theLayer = TECDetId(id).wheel();  // or ring or petal ?
00411           if(TECDetId(id).stereo()==1) { isStereoDet = true ; }
00412         } else {
00413           LogDebug("") << " UHOH BIG PROBLEM - Unrecognized SI Layer" ;
00414           LogDebug("") << " Det "<< theDet << " Lay " << theLayer ;
00415           assert(1!=1) ;
00416         }
00417 
00418         if ((endcap  &&  stereo && (theDet=="TID" || theDet== "TEC") && isStereoDet ) ||
00419             (endcap  &&  !stereo && (theDet=="TID" || theDet== "TEC") && !isStereoDet )  ||
00420             (!endcap  && stereo && (theDet=="TIB" || theDet=="TOB") &&  isStereoDet )    ||
00421             (!endcap  &&  !stereo && (theDet=="TIB" || theDet=="TOB" )&& !isStereoDet )
00422             ) {  
00423               
00424 
00425           hitPointersOut.push_back(&(*hit));
00426 
00427         } // end if this is the right subdetector
00428       } // end loop over hits
00429     } // end if this detector id doesn't have too many hits on it
00430   } // end loop over detector ids
00431 }
00432 
00433 // select all matched hits for now
00434 
00435 void SiStripElectronAlgo::coarseMatchedHitSelection(std::vector<const SiStripMatchedRecHit2D*>& coarseMatchedHitPointersOut)
00436 {
00437   
00438   // Loop over the detector ids
00439   SiStripMatchedRecHit2DCollection::const_iterator itdet = matchedHits_p_->begin(), eddet = matchedHits_p_->end();
00440   for (; itdet != eddet; ++itdet) {
00441     
00442     // Get the hits on this detector id
00443     SiStripMatchedRecHit2DCollection::DetSet hits = *itdet ;
00444     
00445     // Count the number of hits on this detector id
00446     unsigned int numberOfHits = 0;
00447     for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = hits.begin();  hit != hits.end();  ++hit) {
00448       if ( !((hit->geographicalId()).subdetId() == StripSubdetector::TIB) &&
00449            !( (hit->geographicalId()).subdetId() == StripSubdetector::TOB )) { break;}
00450       numberOfHits++;
00451       if (numberOfHits > maxHitsOnDetId_) { break; }
00452     }
00453     
00454     // Only take the hits if there aren't too many
00455     if (numberOfHits <= maxHitsOnDetId_) {
00456       for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = hits.begin();  hit != hits.end();  ++hit) {
00457         if(!(*hit).isValid()) {
00458           LogDebug("") << " InValid hit skipped in coarseMatchedHitSelection " << std::endl ;
00459           continue ;
00460         }
00461         if ( !((hit->geographicalId()).subdetId() == StripSubdetector::TIB) &&
00462              !( (hit->geographicalId()).subdetId() == StripSubdetector::TOB )) { break;}
00463         
00464         coarseMatchedHitPointersOut.push_back(&(*hit));
00465       } // end loop over hits
00466       
00467     } // end if this detector id doesn't have too many hits on it
00468   } // end loop over detector ids
00469   
00470 
00471 }// end of matchedHitSelection
00472 
00473 
00474 
00475 
00476 // projects a phi band of width phiBandWidth_ from supercluster into tracker (given a chargeHypothesis)
00477 // fills *_pos_ or *_neg_ member data with the results
00478 // returns true iff the electron/positron passes cuts
00479 bool SiStripElectronAlgo::projectPhiBand(float chargeHypothesis, const reco::SuperClusterRef& superclusterIn)
00480 {
00481   // This algorithm projects a phi band into the tracker three times:
00482   // (a) for all stereo hits, (b) for barrel rphi hits, and (c) for
00483   // endcap zphi hits.  While accumulating stereo hits in step (a),
00484   // we fit r vs z to a line.  This resolves the ambiguity in z for
00485   // rphi hits and the ambiguity in r for zphi hits.  We can then cut
00486   // on the z of rphi hits (a little wider than one strip length),
00487   // and we can convert the z of zphi hits into r to apply the phi
00488   // band cut.  (We don't take advantage of the endcap strips'
00489   // segmentation in r.)
00490   // 
00491   // As we project a phi band into the tracker, we count hits within
00492   // that band and performs a linear fit for phi vs r.  The number of
00493   // hits and reduced chi^2 from the fit are used to select a good
00494   // candidate.
00495 
00496   // set limits on Si hits - smallest error that makes sense = 1 micron ?
00497   static const double rphiHitSmallestError = 0.0001 ;
00498   static const double stereoHitSmallestError = 0.0001 ;
00499   //not used since endcap is not implemented  
00500   // static const double zphiHitSmallestError = 0.0001 ;
00501 
00502 
00503   static const double stereoPitchAngle = 0.100 ; // stereo angle in rad
00504   static const double cos2SiPitchAngle = cos(stereoPitchAngle)*cos(stereoPitchAngle) ;
00505   static const double sin2SiPitchAngle = sin(stereoPitchAngle)*sin(stereoPitchAngle) ;
00506   // overall misalignment fudge to be added in quad to position errors.
00507   // this is a rough approx to values reported in tracking meet 5/16/2007
00508   static const double rphiErrFudge = 0.0200 ;
00509   static const double stereoErrFudge = 0.0200 ;
00510 
00511   // max chi2 of a hit on an SiDet relative to the prediction
00512   static const double chi2HitMax = 25.0 ;
00513 
00514   // Minimum number of hits to consider on a candidate
00515   static const unsigned int nHitsLeftMinimum = 3  ;
00516 
00517   // Create and fill vectors of pointers to hits
00518   std::vector<const SiStripRecHit2D*> stereoHits;
00519   std::vector<const SiStripRecHit2D*> rphiBarrelHits;
00520   std::vector<const SiStripRecHit2D*> zphiEndcapHits;
00521 
00522   //                                 stereo? endcap?
00523   coarseHitSelection(stereoHits,     true,   false);
00524 
00525   // skip endcap stereo for now
00526   //  LogDebug("") << " Getting endcap stereo hits " ;
00527   // coarseHitSelection(stereoHits,     true,   true);
00528   std::ostringstream debugstr1;
00529   debugstr1 << " Getting barrel rphi hits " << " \n" ;
00530 
00531   coarseHitSelection(rphiBarrelHits, false,  false);
00532 
00533   //  LogDebug("") << " Getting endcap zphi hits " ;
00534   //  coarseHitSelection(zphiEndcapHits, false,  true);
00535 
00536   debugstr1 << " Getting matched hits "  << " \n" ;
00537   std::vector<const SiStripMatchedRecHit2D*> matchedHits;
00538   coarseMatchedHitSelection(matchedHits);
00539 
00540 
00541   // Determine how to project from the supercluster into the tracker
00542   double energy = superclusterIn->energy();
00543   double pT = energy * superclusterIn->position().rho()/sqrt(superclusterIn->x()*superclusterIn->x() +
00544                                                              superclusterIn->y()*superclusterIn->y() +
00545                                                              superclusterIn->z()*superclusterIn->z());
00546   // cf Jackson p. 581-2, a little geometry
00547   double phiVsRSlope = -3.00e-3 * chargeHypothesis * magneticField_p_->inTesla(GlobalPoint(superclusterIn->x(), superclusterIn->y(), 0.)).z() / pT / 2.;
00548 
00549   // Shorthand for supercluster radius, z
00550   const double scr = superclusterIn->position().rho();
00551   const double scz = superclusterIn->position().z();
00552 
00553   // These are used to fit all hits to a line in phi(r)
00554   std::vector<bool> uselist;
00555   std::vector<double> rlist, philist, w2list;
00556   std::vector<int> typelist;  // stereo = 0, rphi barrel = 1, and zphi disk = 2 (only used in this function)
00557   std::vector<const SiStripRecHit2D*> hitlist;
00558 
00559   std::vector<bool> matcheduselist;
00560   std::vector<const SiStripMatchedRecHit2D*> matchedhitlist;
00561 
00562   // These are used to fit the stereo hits to a line in z(r), constrained to pass through supercluster
00563   double zSlopeFitNumer = 0.;
00564   double zSlopeFitDenom = 0.;
00565 
00566   
00567   debugstr1 << " There are a total of " << stereoHits.size()  << " stereoHits in this event " << " \n" 
00568             << " There are a total of " <<  rphiBarrelHits.size() << " rphiBarrelHits in this event " << " \n"
00569             << " There are a total of " <<  zphiEndcapHits.size() << " zphiEndcapHits in this event " << " \n\n";
00570 
00571 
00572   LogDebug("") << debugstr1.str() ;
00573   
00575 
00576 
00577   // Loop over all matched hits
00578   // make a list of good matched rechits.
00579   // in the stereo and rphi loops check to see if the hit is associated with a matchedhit
00580   LogDebug("") << " Loop over matched hits " << " \n";
00581 
00582   for (std::vector<const SiStripMatchedRecHit2D*>::const_iterator hit = matchedHits.begin() ;
00583        hit != matchedHits.end() ; ++ hit) {
00584     
00585     assert(matchedHitUsed_.find(*hit) != matchedHitUsed_.end());
00586 
00587     if (!matchedHitUsed_[*hit]) {
00588 
00589       // Calculate the 3-D position of this hit
00590       GlobalPoint position = tracker_p_->idToDet((*hit)->geographicalId())->surface().toGlobal((*hit)->localPosition());
00591       double r = sqrt(position.x()*position.x() + position.y()*position.y());
00592       double phi = unwrapPhi(position.phi() - superclusterIn->position().phi());  // phi is relative to supercluster
00593       double z = position.z();
00594 
00595       // Cut a triangle in the z-r plane using the supercluster's eta/dip angle information
00596       // and the fact that the electron originated *near* the origin
00597       if ((-originUncertainty_ + (scz + originUncertainty_)*(r/scr)) < z  &&  z < (originUncertainty_ + (scz - originUncertainty_)*(r/scr))) {
00598          
00599         // Cut a narrow band around the supercluster's projection in phi
00600         if (unwrapPhi((r-scr)*phiVsRSlope - phiBandWidth_) < phi  &&  phi < unwrapPhi((r-scr)*phiVsRSlope + phiBandWidth_)) {
00601          
00602 
00603           matcheduselist.push_back(true);
00604           matchedhitlist.push_back(*hit);
00605           
00606 
00607 
00608           // Use this hit to fit z(r)
00609           zSlopeFitNumer += -(scr - r) * (z - scz);
00610           zSlopeFitDenom += (scr - r) * (scr - r);
00611             
00612         } // end cut on phi band
00613       } // end cut on electron originating *near* the origin
00614     } // end assign disjoint sets of hits to electrons
00615   } // end loop over matched hits
00616 
00617   // Calculate the linear fit for z(r)
00618   double zVsRSlope;
00619   if (zSlopeFitDenom > 0.) {
00620     zVsRSlope = zSlopeFitNumer / zSlopeFitDenom;
00621   }
00622   else {
00623     // zVsRSlope assumes electron is from origin if there were no stereo hits
00624     zVsRSlope = scz/scr;
00625   }
00626 
00627   //  // Loop over all stereo hits
00628   LogDebug("") << " Loop over stereo hits" << " \n";
00629 
00630   // check if the stereo hit is matched to one of the matched hit
00631   unsigned int numberOfStereoHits = 0;
00632   for (std::vector<const SiStripRecHit2D*>::const_iterator hit = stereoHits.begin();  hit != stereoHits.end();  ++hit) {
00633     assert(hitUsed_.find(*hit) != hitUsed_.end());
00634 
00635     // Calculate the 3-D position of this hit
00636     GlobalPoint position = tracker_p_->idToDet((*hit)->geographicalId())->surface().toGlobal((*hit)->localPosition());
00637     double r_stereo = sqrt(position.x()*position.x() + position.y()*position.y());
00638     double phi_stereo = unwrapPhi(position.phi() - superclusterIn->position().phi());  // phi is relative to supercluster
00639     //    double z_stereo = position.z();
00640 
00641     // stereo is has a pitch of 100 mrad - so consider both components
00642     double r_stereo_err = sqrt((*hit)->localPositionError().xx()*cos2SiPitchAngle +
00643                                (*hit)->localPositionError().yy()*sin2SiPitchAngle ) ; 
00644 
00645 
00646 
00647     // make sure that the error for this hit is sensible, ie > 1 micron
00648     // otherwise skip this hit
00649     if(r_stereo_err > stereoHitSmallestError ) {
00650       r_stereo_err = sqrt(r_stereo_err*r_stereo_err+stereoErrFudge*stereoErrFudge);
00651  
00652       const SiStripRecHit2D::ClusterRef & stereocluster=(*hit)->cluster();
00653     
00654       bool thisHitIsMatched = false ;
00655 
00656       if (!hitUsed_[*hit]) {
00657  
00658         unsigned int matcheduselist_size = matcheduselist.size();
00659         for (unsigned int i = 0;  i < matcheduselist_size;  i++) {
00660           if (matcheduselist[i]) {
00661             const SiStripRecHit2D::ClusterRef &  mystereocluster = matchedhitlist[i]->stereoHit()->cluster();
00662             if( stereocluster == mystereocluster ) {
00663               thisHitIsMatched = true ;
00664               //    LogDebug("")<< "     This hit is matched " << tracker_p_->idToDet(matchedhitlist[i]->stereoHit()->geographicalId())->surface().toGlobal(matchedhitlist[i]->stereoHit()->localPosition()) << std::endl;
00665               //      break ;
00666             }
00667           } // check if matcheduselist okay 
00668         }// loop over matched hits 
00669       
00670         if(thisHitIsMatched) {
00671           // Use this hit to fit phi(r)
00672           uselist.push_back(true);
00673           rlist.push_back(r_stereo);
00674           philist.push_back(phi_stereo);
00675           w2list.push_back(1./(r_stereo_err/r_stereo)/(r_stereo_err/r_stereo));  // weight**2 == 1./uncertainty**2
00676           typelist.push_back(0);
00677           hitlist.push_back(*hit);
00678         } // thisHitIsMatched
00679       } //  if(!hitUsed)
00680 
00681     } //  end of check on hit position error size 
00682 
00683   } // end loop over stereo hits
00684   
00685   LogDebug("") << " There are " << uselist.size()  << " good hits after stereo loop "  ;
00686  
00687 
00688   // Loop over barrel rphi hits
00689   LogDebug("") << " Looping over barrel rphi hits " ;
00690   unsigned int rphiMatchedCounter = 0 ;
00691   unsigned int rphiUnMatchedCounter = 0 ;
00692   unsigned int numberOfBarrelRphiHits = 0;
00693   for (std::vector<const SiStripRecHit2D*>::const_iterator hit = rphiBarrelHits.begin();  hit != rphiBarrelHits.end();  ++hit) {
00694     assert(hitUsed_.find(*hit) != hitUsed_.end());
00695     // Calculate the 2.5-D position of this hit
00696     GlobalPoint position = tracker_p_->idToDet((*hit)->geographicalId())->surface().toGlobal((*hit)->localPosition());
00697     double r = sqrt(position.x()*position.x() + position.y()*position.y());
00698     double phi = unwrapPhi(position.phi() - superclusterIn->position().phi());  // phi is relative to supercluster
00699     double z = position.z();
00700     double r_mono_err = sqrt((*hit)->localPositionError().xx()) ;
00701 
00702     // only consider hits with errors that make sense
00703     if( r_mono_err > rphiHitSmallestError) {
00704       // inflate the reported error
00705       r_mono_err=sqrt(r_mono_err*r_mono_err+rphiErrFudge*rphiErrFudge);
00706 
00707       const SiStripRecHit2D::ClusterRef &  monocluster=(*hit)->cluster();
00708     
00709 
00710       if (!hitUsed_[*hit]) {
00711         if( (tracker_p_->idToDetUnit((*hit)->geographicalId())->type().subDetector() == GeomDetEnumerators::TIB  && 
00712              (TIBDetId((*hit)->geographicalId()).layer()==1 || TIBDetId((*hit)->geographicalId()).layer()==2)) || 
00713             (tracker_p_->idToDetUnit((*hit)->geographicalId())->type().subDetector() == GeomDetEnumerators::TOB  
00714              && (TOBDetId((*hit)->geographicalId()).layer()==1 || TOBDetId((*hit)->geographicalId()).layer()==2)) ) {
00715           bool thisHitIsMatched = false ;
00716           unsigned int matcheduselist_size = matcheduselist.size();
00717           for (unsigned int i = 0;  i < matcheduselist_size;  i++) {
00718             if (matcheduselist[i]) {
00719               const SiStripRecHit2D::ClusterRef &  mymonocluster = matchedhitlist[i]->monoHit()->cluster();
00720               if( monocluster == mymonocluster ) {
00721                 thisHitIsMatched = true ;
00722               } 
00723             } // check if matcheduselist okay 
00724           }// loop over matched hits 
00725         
00726         
00727           if( thisHitIsMatched ) {
00728             // Use this hit to fit phi(r)
00729             uselist.push_back(true);
00730             rlist.push_back(r);
00731             philist.push_back(phi);
00732             w2list.push_back(1./(r_mono_err/r)/(r_mono_err/r));  // weight**2 == 1./uncertainty**2
00733             typelist.push_back(1);
00734             hitlist.push_back(*hit);
00735             rphiMatchedCounter++;
00736           } // end of matched hit check
00737 
00738         } else {
00739 
00740 
00741           // The expected z position of this hit, according to the z(r) fit
00742           double zFit = zVsRSlope * (r - scr) + scz;
00743         
00744           // Cut on the Z of the strip
00745           // TIB strips are 11 cm long, TOB strips are 19 cm long (can I get these from a function?)
00746           if ((tracker_p_->idToDetUnit((*hit)->geographicalId())->type().subDetector() == GeomDetEnumerators::TIB  && 
00747                std::abs(z - zFit) < 12.)  ||
00748               (tracker_p_->idToDetUnit((*hit)->geographicalId())->type().subDetector() == GeomDetEnumerators::TOB  && 
00749                std::abs(z - zFit) < 20.)    ) {
00750           
00751             // Cut a narrow band around the supercluster's projection in phi
00752             if (unwrapPhi((r-scr)*phiVsRSlope - phiBandWidth_) < phi  &&  phi < unwrapPhi((r-scr)*phiVsRSlope + phiBandWidth_)) {
00753             
00754               // Use this hit to fit phi(r)
00755               uselist.push_back(true);
00756               rlist.push_back(r);
00757               philist.push_back(phi);
00758               w2list.push_back(1./(r_mono_err/r)/(r_mono_err/r));  // weight**2 == 1./uncertainty**2
00759               typelist.push_back(1);
00760               hitlist.push_back(*hit);
00761               rphiUnMatchedCounter++;
00762             
00763             } // end cut on phi band
00764           } // end cut on strip z
00765         } // loop over TIB/TOB layer 1,2 
00766       } // end assign disjoint sets of hits to electrons
00767     } // end of check on rphi hit position error size
00768   } // end loop over barrel rphi hits
00769 
00770   LogDebug("") << " There are " << rphiMatchedCounter <<" matched rphi hits"; 
00771   LogDebug("") << " There are " << rphiUnMatchedCounter <<" unmatched rphi hits";
00772   LogDebug("") << " There are " << uselist.size() << " good stereo+rphi hits " ;
00773 
00774 
00775 
00776 
00777 
00778 
00780 
00781   // Loop over endcap zphi hits
00782   LogDebug("") << " Looping over barrel zphi hits " ;
00783 
00784 
00785   unsigned int numberOfEndcapZphiHits = 0;
00786   for (std::vector<const SiStripRecHit2D*>::const_iterator hit = zphiEndcapHits.begin();  
00787        hit != zphiEndcapHits.end();  ++hit) {
00788     assert(hitUsed_.find(*hit) != hitUsed_.end());
00789     if (!hitUsed_[*hit]) {
00790       // Calculate the 2.5-D position of this hit
00791       GlobalPoint position = tracker_p_->idToDet((*hit)->geographicalId())->surface().toGlobal((*hit)->localPosition());
00792       double z = position.z();
00793       double phi = unwrapPhi(position.phi() - superclusterIn->position().phi());  // phi is relative to supercluster
00794       //      double r=sqrt(position.x()*position.x()+position.y()*position.y()) ;
00795 
00796       // The expected r position of this hit, according to the z(r) fit
00797       double rFit = (z - scz)/zVsRSlope + scr;
00798 
00799       // I don't know the r widths of the endcap strips, otherwise I
00800       // would apply a cut on r similar to the rphi z cut
00801 
00802       // Cut a narrow band around the supercluster's projection in phi,
00803       // and be sure the hit isn't in a reflected band (extrapolation of large z differences into negative r)
00804       if (rFit > 0.  &&
00805           unwrapPhi((rFit-scr)*phiVsRSlope - phiBandWidth_) < phi  &&  phi < unwrapPhi((rFit-scr)*phiVsRSlope + phiBandWidth_)) {
00806 
00807         // Use this hit to fit phi(r)
00808         uselist.push_back(true);
00809         rlist.push_back(rFit);
00810         philist.push_back(phi);
00811         w2list.push_back(1./(0.05/rFit)/(0.05/rFit));  // weight**2 == 1./uncertainty**2
00812         typelist.push_back(2);
00813         hitlist.push_back(*hit);
00814 
00815       } // end cut on phi band
00816     } // end assign disjoint sets of hits to electrons
00817   } // end loop over endcap zphi hits
00818  
00819   LogDebug("") << " There are " << uselist.size() << " good stereo+rphi+zphi hits " ;
00821 
00822   std::ostringstream debugstr5;
00823   debugstr5 << " List of hits before uniqification " << " \n" ;
00824   for (unsigned int i = 0;  i < uselist.size();  i++) {
00825     if ( uselist[i] ) {
00826       const SiStripRecHit2D* hit = hitlist[i];
00827       debugstr5 << " DetID " << ((hit)->geographicalId()).rawId()
00828                 << " R " << rlist[i] 
00829                 << " Phi " << philist[i]
00830                 << " Weight " << w2list[i] 
00831                 << " PhiPred " << (rlist[i]-scr)*phiVsRSlope  
00832                 << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
00833                 << " \n" ;
00834     }
00835   }
00836   debugstr5 << " \n\n\n" ;
00837   
00838   debugstr5 << " Counting number of unique detectors " << " \n" ;
00839 
00840   debugstr5 << " These are all the detectors with hits " << " \n";
00841   // Loop over hits, and find the best hit for each SiDetUnit - keep only those
00842   // get a list of DetIds in uselist
00843   std::vector<uint32_t> detIdList ;
00844 
00845   for (unsigned int i = 0;  i < uselist.size();  i++) {
00846     if (uselist[i]) {
00847       const SiStripRecHit2D* hit = hitlist[i];
00848       uint32_t detIDUsed = ((hit)->geographicalId()).rawId() ;
00849       debugstr5<< " DetId " << detIDUsed << "\n";
00850       detIdList.push_back(detIDUsed) ;
00851     }
00852   }
00853   debugstr5 << " There are " << detIdList.size() << " hits on detectors \n";
00854   // now sort and then uniq this list of detId
00855   std::sort(detIdList.begin(), detIdList.end());
00856   detIdList.erase(
00857                   std::unique(detIdList.begin(), detIdList.end()), detIdList.end());
00858 
00859   debugstr5 << " There are " << detIdList.size() << " unique detectors \n";
00860 
00861   //now we have a list of unique detectors used
00862 
00863 
00864 
00865   debugstr5 << " Looping over detectors to choose best hit " << " \n";
00866   // Loop over detectors ID and hits to create list of hits on same DetId
00867   for (unsigned int idet = 0 ; idet < detIdList.size() ; idet++ ) {
00868     for (unsigned int i = 0;  i+1 < uselist.size();  i++) {
00869       if (uselist[i]) {
00870         // Get Chi2 of this hit relative to predicted hit
00871         const SiStripRecHit2D* hit1 = hitlist[i];
00872         double r_hit1 = rlist[i];
00873         double phi_hit1 = philist[i];
00874         double w2_hit1 = w2list[i] ;
00875         double phi_pred1 = (r_hit1-scr)*phiVsRSlope ; 
00876         double chi1 = (phi_hit1-phi_pred1)*(phi_hit1-phi_pred1)*w2_hit1;
00877         if(detIdList[idet]== ((hit1)->geographicalId()).rawId()) {
00878           for (unsigned int j = i+1;  j < uselist.size();  j++) {
00879             if (uselist[j]) {
00880               const SiStripRecHit2D* hit2 = hitlist[j];
00881               if(detIdList[idet]== ((hit2)->geographicalId()).rawId()) {
00882                 debugstr5 << " Found 2 hits on same Si Detector " 
00883                           << ((hit2)->geographicalId()).rawId() << "\n";
00884                 // Get Chi2 of this hit relative to predicted hit
00885                 double r_hit2 = rlist[j];
00886                 double phi_hit2 = philist[j];
00887                 double w2_hit2 = w2list[j] ;
00888                 double phi_pred2 = (r_hit2-scr)*phiVsRSlope ; 
00889                 double chi2 = (phi_hit2-phi_pred2)*(phi_hit2-phi_pred2)*w2_hit2;
00890                 debugstr5 << " Chi1 " << chi1 << " Chi2 " << chi2 <<"\n";
00891                 if( chi1< chi2 ){
00892                   uselist[j] = 0;
00893                 }else{
00894                   uselist[i] = 0;
00895                 }
00896 
00897               } // end of Det check
00898             } // end of j- usehit check
00899           } // end of j hit list loop
00900           
00901         } // end of detector check
00902       } // end of uselist check
00903     } // end of i-hit list loop
00904 
00905 
00906   } // end of DetId loop
00907 
00908 
00909   
00910   // now let's through out hits with a predicted chi > chi2HitMax 
00911   for ( unsigned int i = 0;  i < uselist.size();  i++ ) { 
00912     if ( uselist[i] ) { 
00913       const SiStripRecHit2D* hit = hitlist[i];
00914       double localchi2 = (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i] ;
00915       if(localchi2 > chi2HitMax ) {
00916         debugstr5 << " Throwing out "
00917                   <<" DetID " << ((hit)->geographicalId()).rawId()
00918                   << " R " << rlist[i] 
00919                   << " Phi " << philist[i]
00920                   << " Weight " << w2list[i] 
00921                   << " PhiPred " << (rlist[i]-scr)*phiVsRSlope  
00922                   << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
00923                   << " \n" ;
00924         uselist[i]=0 ;
00925       }
00926     }
00927   }
00928 
00929   debugstr5 << " List of hits after uniqification " << " \n" ;
00930   for (unsigned int i = 0;  i < uselist.size();  i++) {
00931     if ( uselist[i] ) {
00932       const SiStripRecHit2D* hit = hitlist[i];
00933       debugstr5 << " DetID " << ((hit)->geographicalId()).rawId()
00934                 << " R " << rlist[i] 
00935                 << " Phi " << philist[i]
00936                 << " Weight " << w2list[i] 
00937                 << " PhiPred " << (rlist[i]-scr)*phiVsRSlope  
00938                 << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
00939                 << " \n" ;
00940     }
00941   }
00942   debugstr5 << " \n\n\n" ;
00943 
00944 
00945   // need to check that there are more than 2 hits left here!
00946   unsigned int nHitsLeft =0;
00947   for (unsigned int i = 0;  i < uselist.size();  i++) {
00948     if ( uselist[i] ) {
00949       nHitsLeft++;
00950     }
00951   }
00952   if(nHitsLeft < nHitsLeftMinimum ) {
00953     debugstr5 << " Too few hits "<< nHitsLeft 
00954               << " left - return false " << " \n";
00955     return false ;
00956   }
00957   LogDebug("") << debugstr5.str();
00958 
00960   
00961   // Calculate a linear phi(r) fit and drop hits until the biggest contributor to chi^2 is less than maxNormResid_
00962   bool done = false;
00963   double intercept = 0 ;
00964   double slope = 0 ;
00965   double chi2 = 0;
00966 
00967   std::ostringstream debugstr4;
00968   debugstr4 <<" Calc of phi(r) "<<" \n";
00969 
00970   while (!done) {
00971     // Use an iterative update of <psi>, <r> etc instead in an
00972     // attempt to minize divisions of  large numbers
00973     // The linear fit  psi = slope*r + intercept
00974     //             slope is (<psi*r>-<psi>*<r>) / (<r^2>-<r>^2>)
00975     //             intercept is <psi> - slope*<r>
00976     
00977     double phiBar   = 0.;
00978     double phiBarOld   = 0.;
00979     double rBar  = 0.;
00980     double rBarOld  = 0.;
00981     double r2Bar = 0.;
00982     double r2BarOld = 0.;
00983     double phirBar = 0.;
00984     double phirBarOld = 0.;
00985     double totalWeight = 0.;
00986     double totalWeightOld = 0.; 
00987     unsigned int uselist_size = uselist.size();
00988     for (unsigned int i = 0;  i < uselist_size;  i++) {
00989       if (uselist[i]) {
00990         double r = rlist[i];
00991         double phi = philist[i];
00992         double weight2 = w2list[i];
00993         
00994         totalWeight = totalWeightOld + weight2 ;
00995         
00996         //weight2 is 1/sigma^2. Typically sigma is 100micron pitch
00997         // over root(12) = 30 microns so weight2 is about 10^5-10^6
00998         
00999         double totalWeightRatio = totalWeightOld/totalWeight ;
01000         double localWeightRatio = weight2/totalWeight ;
01001         
01002         phiBar = phiBarOld*totalWeightRatio + phi*localWeightRatio ; 
01003         rBar = rBarOld*totalWeightRatio + r*localWeightRatio ; 
01004         r2Bar= r2BarOld*totalWeightRatio + r*r*localWeightRatio ; 
01005         phirBar = phirBarOld*totalWeightRatio + phi*r*localWeightRatio ;
01006         
01007         totalWeightOld = totalWeight ;
01008         phiBarOld = phiBar ;
01009         rBarOld = rBar ;
01010         r2BarOld = r2Bar ;
01011         phirBarOld = phirBar ;  
01012         
01013         debugstr4 << " totalWeight " << totalWeight 
01014                   << " totalWeightRatio " << totalWeightRatio
01015                   << " localWeightRatio "<< localWeightRatio
01016                   << " phiBar " << phiBar
01017                   << " rBar " << rBar 
01018                   << " r2Bar " << r2Bar
01019                   << " phirBar " << phirBar 
01020                   << " \n ";
01021 
01022 
01023       } // end of use hit loop
01024     } // end of hit loop to calculate a linear fit
01025     slope = (phirBar-phiBar*rBar)/(r2Bar-rBar*rBar);
01026     intercept = phiBar-slope*rBar ;
01027 
01028     debugstr4 << " end of loop slope " << slope 
01029               << " intercept " << intercept << " \n";
01030 
01031     // Calculate chi^2
01032     chi2 = 0.;
01033     double biggest_normresid = -1.;
01034     unsigned int biggest_index = 0;
01035     for (unsigned int i = 0;  i < uselist_size;  i++) {
01036       if (uselist[i]) {
01037         double r = rlist[i];
01038         double phi = philist[i];
01039         double weight2 = w2list[i];
01040 
01041         double normresid = weight2 * (phi - intercept - slope*r)*(phi - intercept - slope*r);
01042         chi2 += normresid;
01043 
01044         if (normresid > biggest_normresid) {
01045           biggest_normresid = normresid;
01046           biggest_index = i;
01047         }
01048       }
01049     } // end loop over hits to calculate chi^2 and find its biggest contributer
01050 
01051     if (biggest_normresid > maxNormResid_) {
01052       debugstr4 << "Dropping hit from fit due to Chi2 " << " \n" ;
01053       const SiStripRecHit2D* hit = hitlist[biggest_index];
01054       debugstr4 << " DetID " << ((hit)->geographicalId()).rawId()
01055                 << " R " << rlist[biggest_index]
01056                 << " Phi " << philist[biggest_index]
01057                 << " Weight " << w2list[biggest_index]
01058                 << " PhiPred " << (rlist[biggest_index]-scr)*phiVsRSlope  
01059                 << " Chi2 " << (philist[biggest_index]-(rlist[biggest_index]-scr)*phiVsRSlope)*(philist[biggest_index]-(rlist[biggest_index]-scr)*phiVsRSlope)*w2list[biggest_index] 
01060                 << " normresid " <<  biggest_normresid
01061                 << " \n\n";
01062       uselist[biggest_index] = false;
01063     }
01064     else {
01065       done = true;
01066     }
01067   } // end loop over trial fits
01068 
01069   debugstr4 <<" Fit " 
01070             << " Intercept  " << intercept
01071             << " slope " << slope 
01072             << " chi2 " << chi2
01073             << " \n" ;
01074 
01075   LogDebug("") <<   debugstr4.str() ;
01076 
01077 
01078   // Now we have intercept, slope, and chi2; uselist to tell us which hits are used, and hitlist for the hits
01079 
01080   // Identify the innermost hit
01081   const SiStripRecHit2D* innerhit = (SiStripRecHit2D*)(0);
01082   double innerhitRadius = -1.;  // meaningless until innerhit is defined
01083 
01084   // Copy hits into an OwnVector, which we put in the TrackCandidate
01085   std::vector<const TrackingRecHit*> outputHits;
01086   // Reference rphi and stereo hits into RefVectors, which we put in the SiStripElectron
01087   std::vector<SiStripRecHit2D> outputRphiHits;
01088   std::vector<SiStripRecHit2D> outputStereoHits;
01089 
01090   typedef edm::Ref<SiStripRecHit2DCollection,SiStripRecHit2D> SiStripRecHit2DRef;
01091 
01092 
01093   for (unsigned int i = 0;  i < uselist.size();  i++) {
01094     if (uselist[i]) {
01095       double r = rlist[i];
01096       const SiStripRecHit2D* hit = hitlist[i];
01097 
01098       // Keep track of the innermost hit
01099       if (innerhit == (SiStripRecHit2D*)(0)  ||  r < innerhitRadius) {
01100         innerhit = hit;
01101         innerhitRadius = r;
01102       }
01103          
01104       if (typelist[i] == 0) {
01105         numberOfStereoHits++;
01106 
01107         // Copy this hit for the TrajectorySeed
01108         outputHits.push_back(hit);
01109         outputStereoHits.push_back(*hit);
01110       }
01111       else if (typelist[i] == 1) {
01112         numberOfBarrelRphiHits++;
01113 
01114         // Copy this hit for the TrajectorySeed
01115         outputHits.push_back(hit);
01116         outputRphiHits.push_back(*hit);
01117       }
01118       else if (typelist[i] == 2) {
01119         numberOfEndcapZphiHits++;
01120 
01121         // Copy this hit for the TrajectorySeed
01122         outputHits.push_back(hit);
01123         outputRphiHits.push_back(*hit);
01124       }
01125     }
01126   } // end loop over all hits, after having culled the ones with big residuals
01127 
01128   unsigned int totalNumberOfHits = numberOfStereoHits + numberOfBarrelRphiHits + numberOfEndcapZphiHits;
01129   double reducedChi2 = (totalNumberOfHits > 2 ? chi2 / (totalNumberOfHits - 2) : 1e10);
01130 
01131   // Select this candidate if it passes minHits_ and maxReducedChi2_ cuts
01132   if (totalNumberOfHits >= minHits_  &&  reducedChi2 <= maxReducedChi2_) {
01133     // GlobalTrajectoryParameters evaluated at the position of the innerhit
01134     GlobalPoint position = tracker_p_->idToDet(innerhit->geographicalId())->surface().toGlobal(innerhit->localPosition());
01135 
01136     // Use our phi(r) linear fit to correct pT (pT is inversely proportional to slope)
01137     // (By applying a correction instead of going back to first
01138     // principles, any correction to the phiVsRSlope definition
01139     // above will be automatically propagated here.)
01140     double correct_pT = pT * phiVsRSlope / slope;
01141 
01142     // Our phi(r) linear fit returns phi relative to the supercluster phi for a given radius
01143     double phi = intercept + slope*sqrt(position.x()*position.x() + position.y()*position.y()) + superclusterIn->position().phi();
01144 
01145     // Our z(r) linear fit gives us a better measurement of eta/dip angle
01146     double pZ = correct_pT * zVsRSlope;
01147 
01148     // Now we can construct a trajectory momentum out of linear fits to hits
01149     GlobalVector momentum = GlobalVector(correct_pT*cos(phi), correct_pT*sin(phi), pZ);
01150 
01151     if (chargeHypothesis > 0.) {
01152       redchi2_pos_ = chi2 / double(totalNumberOfHits - 2);
01153       position_pos_ = position;
01154       momentum_pos_ = momentum;
01155       innerhit_pos_ = innerhit;
01156       outputHits_pos_ = outputHits;
01157       outputRphiHits_pos_ = outputRphiHits;
01158       outputStereoHits_pos_ = outputStereoHits;
01159       phiVsRSlope_pos_ = phiVsRSlope;
01160       slope_pos_ = slope;
01161       intercept_pos_ = intercept;
01162       chi2_pos_ = chi2;
01163       ndof_pos_ = totalNumberOfHits - 2;
01164       correct_pT_pos_ = correct_pT;
01165       pZ_pos_ = pZ;
01166       zVsRSlope_pos_ = zVsRSlope;
01167       numberOfStereoHits_pos_ = numberOfStereoHits;
01168       numberOfBarrelRphiHits_pos_ = numberOfBarrelRphiHits;
01169       numberOfEndcapZphiHits_pos_ = numberOfEndcapZphiHits;
01170     }
01171     else {
01172       redchi2_neg_ = chi2 / double(totalNumberOfHits - 2);
01173       position_neg_ = position;
01174       momentum_neg_ = momentum;
01175       innerhit_neg_ = innerhit;
01176       outputHits_neg_ = outputHits;
01177       outputRphiHits_neg_ = outputRphiHits;
01178       outputStereoHits_neg_ = outputStereoHits;
01179       phiVsRSlope_neg_ = phiVsRSlope;
01180       slope_neg_ = slope;
01181       intercept_neg_ = intercept;
01182       chi2_neg_ = chi2;
01183       ndof_neg_ = totalNumberOfHits - 2;
01184       correct_pT_neg_ = correct_pT;
01185       pZ_neg_ = pZ;
01186       zVsRSlope_neg_ = zVsRSlope;
01187       numberOfStereoHits_neg_ = numberOfStereoHits;
01188       numberOfBarrelRphiHits_neg_ = numberOfBarrelRphiHits;
01189       numberOfEndcapZphiHits_neg_ = numberOfEndcapZphiHits;
01190     }
01191 
01192     LogDebug("") << " return from projectFindBand with True \n" << std::endl ;
01193     return true;
01194   } // end if this is a good electron candidate
01195 
01196   // Signal for a failed electron candidate
01197   LogDebug("") << " return from projectFindBand with False \n" << std::endl ;
01198   return false;
01199 }
01200 
01201 //
01202 // const member functions
01203 //
01204 
01205 //
01206 // static member functions
01207 //
01208