CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/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.40 2012/01/16 09:28:17 innocent 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 
00022 #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
00023 #include "RecoTracker/TrackProducer/interface/TrackingRecHitLessFromGlobalPosition.h"
00024 
00025 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00026 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00027 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00028 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00029 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00030 
00031 
00032 //
00033 // constants, enums and typedefs
00034 //
00035 
00036 //
00037 // static data member definitions
00038 //
00039 namespace {
00040   struct CompareHits {
00041     CompareHits( const TrackingRecHitLessFromGlobalPosition& iLess ) :
00042       less_(iLess) {}
00043     bool operator()(const TrackingRecHit* iLHS, const TrackingRecHit* iRHS) {
00044       return less_(*iLHS,*iRHS);
00045     }
00046          
00047     TrackingRecHitLessFromGlobalPosition less_;
00048   };
00049 }
00050 
00051 
00052 //
00053 // constructors and destructor
00054 //
00055 SiStripElectronAlgo::SiStripElectronAlgo(unsigned int maxHitsOnDetId,
00056                                          double originUncertainty,
00057                                          double phiBandWidth,
00058                                          double maxNormResid,
00059                                          unsigned int minHits,
00060                                          double maxReducedChi2)
00061   : maxHitsOnDetId_(maxHitsOnDetId)
00062   , originUncertainty_(originUncertainty)
00063   , phiBandWidth_(phiBandWidth)
00064   , maxNormResid_(maxNormResid)
00065   , minHits_(minHits)
00066   , maxReducedChi2_(maxReducedChi2)
00067 {
00068 }
00069 
00070 // SiStripElectronAlgo::SiStripElectronAlgo(const SiStripElectronAlgo& rhs)
00071 // {
00072 //    // do actual copying here;
00073 // }
00074 
00075 SiStripElectronAlgo::~SiStripElectronAlgo()
00076 {
00077 }
00078 
00079 //
00080 // assignment operators
00081 //
00082 // const SiStripElectronAlgo& SiStripElectronAlgo::operator=(const SiStripElectronAlgo& rhs)
00083 // {
00084 //   //An exception safe implementation is
00085 //   SiStripElectronAlgo temp(rhs);
00086 //   swap(rhs);
00087 //
00088 //   return *this;
00089 // }
00090 
00091 //
00092 // member functions
00093 //
00094 
00095 void SiStripElectronAlgo::prepareEvent(const edm::ESHandle<TrackerGeometry>& tracker,
00096                                        const edm::Handle<SiStripRecHit2DCollection>& rphiHits,
00097                                        const edm::Handle<SiStripRecHit2DCollection>& stereoHits,
00098                                        const edm::Handle<SiStripMatchedRecHit2DCollection>& matchedHits,
00099                                        const edm::ESHandle<MagneticField>& magneticField)
00100 {
00101   LogDebug("") << " In prepareEvent " ; 
00102 
00103   tracker_p_ = tracker.product();
00104   rphiHits_p_ = rphiHits.product();
00105   stereoHits_p_ = stereoHits.product();
00106   matchedHits_p_ = matchedHits.product();
00107   magneticField_p_ = magneticField.product();
00108 
00109   rphiHits_hp_ = &rphiHits;
00110   stereoHits_hp_ = &stereoHits;
00111 
00112   // Keep a table that relates hit pointers to their index (key) in the collections
00113   rphiKey_.clear();
00114   stereoKey_.clear();
00115   // Keep track of which hits have been used already (so a hit is assigned to only one electron)
00116   hitUsed_.clear();
00117   matchedHitUsed_.clear();
00118 
00119   unsigned int counter = 0;
00120   for (SiStripRecHit2DCollection::DataContainer::const_iterator it = rphiHits_p_->data().begin();  it != rphiHits_p_->data().end();  ++it) {
00121     rphiKey_[&(*it)] = counter;
00122     hitUsed_[&(*it)] = false;
00123     counter++;
00124   }
00125 
00126   counter = 0;
00127   for (SiStripRecHit2DCollection::DataContainer::const_iterator it = stereoHits_p_->data().begin();  it != stereoHits_p_->data().end();  ++it) {
00128     stereoKey_[&(*it)] = counter;
00129     hitUsed_[&(*it)] = false;
00130     counter++;
00131   }
00132 
00133   counter = 0;
00134   for (SiStripMatchedRecHit2DCollection::DataContainer::const_iterator it = matchedHits_p_->data().begin();  it != matchedHits_p_->data().end();  ++it) {
00135     matchedKey_[&(*it)] = counter;
00136     matchedHitUsed_[&(*it)] = false;
00137     counter++;
00138   }
00139   
00140   LogDebug("") << " Leaving prepareEvent " ;
00141 }
00142 
00143 
00144 
00145 // returns true iff an electron was found
00146 // inserts electrons and trackcandiates into electronOut and trackCandidateOut
00147 bool SiStripElectronAlgo::findElectron(reco::SiStripElectronCollection& electronOut,
00148                                        TrackCandidateCollection& trackCandidateOut,
00149                                        const reco::SuperClusterRef& superclusterIn)
00150 {
00151   // Try each of the two charge hypotheses, but only take one
00152   bool electronSuccess = projectPhiBand(-1., superclusterIn);
00153   bool positronSuccess = projectPhiBand( 1., superclusterIn);
00154 
00155   // electron hypothesis did better than electron
00156   if ((electronSuccess  &&  !positronSuccess)  ||
00157       (electronSuccess  &&  positronSuccess  &&  redchi2_neg_ <= redchi2_pos_)) {
00158       
00159     // Initial uncertainty for tracking
00160     AlgebraicSymMatrix55 errors;  // makes 5x5 matrix, indexed from (0,0) to (44)
00161     errors(0,0) = 3.;      // uncertainty**2 in 1/momentum
00162     errors(1,1) = 0.01;    // uncertainty**2 in lambda (lambda == pi/2 - polar angle theta)
00163     errors(2,2) = 0.0001;  // uncertainty**2 in phi
00164     errors(3,3) = 0.01;    // uncertainty**2 in x_transverse (where x is in cm)
00165     errors(4,4) = 0.01;    // uncertainty**2 in y_transverse (where y is in cm)
00166 
00167 #ifdef EDM_ML_DEBUG 
00168     // JED Debugging possible double hit sort problem
00169     std::ostringstream debugstr6;
00170     debugstr6 << " HERE BEFORE SORT electron case " << " \n" ;
00171     for (std::vector<const TrackingRecHit*>::iterator itHit=outputHits_neg_.begin(); 
00172          itHit != outputHits_neg_.end(); ++itHit) {
00173       debugstr6 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
00174                 <<"    Local Position: "<<(*itHit)->localPosition()<<" \n"
00175                 <<"    Global Rho:  "
00176                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
00177                 <<" Phi "
00178                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
00179                 << " Z "
00180                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";
00181     }
00182     // end of dump 
00183     
00184     
00185     
00186     // JED call dump 
00187     debugstr6 << " Calling sort alongMomentum " << " \n"; 
00188 #endif   
00189  
00190     std::sort(outputHits_neg_.begin(), outputHits_neg_.end(),
00191               CompareHits(TrackingRecHitLessFromGlobalPosition(((TrackingGeometry*)(tracker_p_)), alongMomentum)));
00192     
00193 #ifdef EDM_ML_DEBUG 
00194     debugstr6 << " Done with sort " << " \n";
00195     
00196     debugstr6 << " HERE AFTER SORT electron case " << " \n";
00197     for (std::vector<const TrackingRecHit*>::iterator itHit=outputHits_neg_.begin(); 
00198          itHit != outputHits_neg_.end(); ++itHit) {
00199       debugstr6 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
00200                 <<"    Local Position: "<<(*itHit)->localPosition()<<" \n"
00201                 <<"    Global Rho:  "
00202                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
00203                 <<" Phi "
00204                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
00205                 << " Z "
00206                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";; 
00207     }
00208     // end of dump 
00209 
00210     LogDebug("")<< debugstr6.str();
00211 #endif
00212     
00213     //create an OwnVector needed by the classes which will be stored to the Event
00214     edm::OwnVector<TrackingRecHit> hits;
00215     for(std::vector<const TrackingRecHit*>::iterator itHit =outputHits_neg_.begin();
00216         itHit != outputHits_neg_.end();
00217         ++itHit) {
00218       hits.push_back( (*itHit)->clone());
00219       if( !(hitUsed_.find(*itHit) != hitUsed_.end()) ) {
00220         LogDebug("") << " Assert failure " ;
00221         assert(hitUsed_.find(*itHit) != hitUsed_.end());
00222       }
00223       hitUsed_[*itHit] = true;
00224     }
00225     
00226     TrajectoryStateOnSurface state(
00227                                    GlobalTrajectoryParameters(position_neg_, momentum_neg_, -1, magneticField_p_),
00228                                    CurvilinearTrajectoryError(errors),
00229                                    tracker_p_->idToDet(innerhit_neg_->geographicalId())->surface());
00230     
00231     
00232     PTrajectoryStateOnDet const & PTraj = trajectoryStateTransform::persistentState(state, innerhit_neg_->geographicalId().rawId());
00233     TrajectorySeed trajectorySeed(PTraj, hits, alongMomentum);
00234     trackCandidateOut.push_back(TrackCandidate(hits, trajectorySeed, PTraj));
00235     
00236     electronOut.push_back(reco::SiStripElectron(superclusterIn,
00237                                                 -1,
00238                                                 outputRphiHits_neg_,
00239                                                 outputStereoHits_neg_,
00240                                                 phiVsRSlope_neg_,
00241                                                 slope_neg_,
00242                                                 intercept_neg_ + superclusterIn->position().phi(),
00243                                                 chi2_neg_,
00244                                                 ndof_neg_,
00245                                                 correct_pT_neg_,
00246                                                 pZ_neg_,
00247                                                 zVsRSlope_neg_,
00248                                                 numberOfStereoHits_neg_,
00249                                                 numberOfBarrelRphiHits_neg_,
00250                                                 numberOfEndcapZphiHits_neg_));
00251     
00252     return true;
00253   }
00254   
00255   // positron hypothesis did better than electron
00256   if ((!electronSuccess  &&  positronSuccess)  ||
00257       (electronSuccess  &&  positronSuccess  &&  redchi2_neg_ > redchi2_pos_)) {
00258       
00259     // Initial uncertainty for tracking
00260     AlgebraicSymMatrix55 errors;  // same as abovr
00261     errors(0,0) = 3.;      // uncertainty**2 in 1/momentum
00262     errors(1,1) = 0.01;    // uncertainty**2 in lambda (lambda == pi/2 - polar angle theta)
00263     errors(2,2) = 0.0001;  // uncertainty**2 in phi
00264     errors(3,3) = 0.01;    // uncertainty**2 in x_transverse (where x is in cm)
00265     errors(4,4) = 0.01;    // uncertainty**2 in y_transverse (where y is in cm)
00266 
00267 #ifdef EDM_ML_DEBUG 
00268     // JED Debugging possible double hit sort problem
00269     std::ostringstream debugstr7;
00270     debugstr7 << " HERE BEFORE SORT Positron case " << " \n";
00271     for (std::vector<const TrackingRecHit*>::iterator itHit = outputHits_pos_.begin(); 
00272          itHit != outputHits_pos_.end(); ++itHit) {
00273       debugstr7 <<" HIT "<<((*itHit)->geographicalId()).rawId()<<" \n"
00274                 <<"    Local Position: "<<(*itHit)->localPosition()<<" \n"
00275                 <<"    Global Rho:  "
00276                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).perp())
00277                 <<" Phi "
00278                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).phi())
00279                 << " Z "
00280                 <<(((TrackingGeometry*)(tracker_p_))->idToDet((*itHit)->geographicalId())->surface().toGlobal((*itHit)->localPosition()).z())<< " \n";
00281     }
00282     // end of dump 
00283     
00284     debugstr7 << " Calling sort alongMomentum " << " \n"; 
00285 #endif
00286     
00287     std::sort(outputHits_pos_.begin(), outputHits_pos_.end(),
00288               CompareHits(TrackingRecHitLessFromGlobalPosition(((TrackingGeometry*)(tracker_p_)), alongMomentum)));
00289     
00290 #ifdef EDM_ML_DEBUG 
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 #endif
00309 
00310     //create an OwnVector needed by the classes which will be stored to the Event
00311     edm::OwnVector<TrackingRecHit> hits;
00312     for(std::vector<const TrackingRecHit*>::iterator itHit =outputHits_pos_.begin();
00313         itHit != outputHits_pos_.end();
00314         ++itHit) {
00315       hits.push_back( (*itHit)->clone());
00316       assert(hitUsed_.find(*itHit) != hitUsed_.end());
00317       hitUsed_[*itHit] = true;
00318     }
00319     
00320     TrajectoryStateOnSurface state(
00321                                    GlobalTrajectoryParameters(position_pos_, momentum_pos_, 1, magneticField_p_),
00322                                    CurvilinearTrajectoryError(errors),
00323                                    tracker_p_->idToDet(innerhit_pos_->geographicalId())->surface());
00324     
00325     
00326     PTrajectoryStateOnDet const & PTraj = trajectoryStateTransform::persistentState(state, innerhit_pos_->geographicalId().rawId());
00327     TrajectorySeed trajectorySeed(PTraj, hits, alongMomentum);
00328     trackCandidateOut.push_back(TrackCandidate(hits, trajectorySeed, PTraj));
00329     
00330     electronOut.push_back(reco::SiStripElectron(superclusterIn,
00331                                                 1,
00332                                                 outputRphiHits_pos_,
00333                                                 outputStereoHits_pos_,
00334                                                 phiVsRSlope_pos_,
00335                                                 slope_pos_,
00336                                                 intercept_pos_ + superclusterIn->position().phi(),
00337                                                 chi2_pos_,
00338                                                 ndof_pos_,
00339                                                 correct_pT_pos_,
00340                                                 pZ_pos_,
00341                                                 zVsRSlope_pos_,
00342                                                 numberOfStereoHits_pos_,
00343                                                 numberOfBarrelRphiHits_pos_,
00344                                                 numberOfEndcapZphiHits_pos_));
00345     
00346     return true;
00347   }
00348   
00349   return false;
00350 }
00351 
00352 // inserts pointers to good hits into hitPointersOut
00353 // selects hits on DetIds that have no more than maxHitsOnDetId_
00354 // selects from stereo if stereo == true, rphi otherwise
00355 // selects from TID or TEC if endcap == true, TIB or TOB otherwise
00356 void SiStripElectronAlgo::coarseHitSelection(std::vector<const SiStripRecHit2D*>& hitPointersOut,
00357                                              bool stereo, bool endcap)
00358 {
00359   // This function is not time-efficient.  If you want to improve the
00360   // speed of the algorithm, you'll probably want to change this
00361   // function.  There may be a more efficienct way to extract hits,
00362   // and it would definitely help to put a geographical cut on the
00363   // DetIds.  (How does one determine the global position of a given
00364   // DetId?  Is tracker_p_->idToDet(id)->surface().toGlobal(LocalPosition(0,0,0))
00365   // expensive?)
00366 
00367   // Loop over the detector ids
00368   SiStripRecHit2DCollection::const_iterator itdet = (stereo ? stereoHits_p_->begin() : rphiHits_p_->begin());
00369   SiStripRecHit2DCollection::const_iterator eddet = (stereo ? stereoHits_p_->end()   : rphiHits_p_->end()  );
00370   for (; itdet != eddet; ++itdet) {
00371     // Get the hits on this detector id
00372     SiStripRecHit2DCollection::DetSet hits = *itdet;
00373     DetId id(hits.detId());
00374 
00375     // Count the number of hits on this detector id
00376     unsigned int numberOfHits = hits.size();
00377       
00378     // Only take the hits if there aren't too many
00379     // (Would it be better to loop only once, fill a temporary list,
00380     // and copy that if numberOfHits <= maxHitsOnDetId_?)
00381     if (numberOfHits <= maxHitsOnDetId_) {
00382       for (SiStripRecHit2DCollection::DetSet::const_iterator hit = hits.begin();  
00383            hit != hits.end();  ++hit) {
00384         // check that hit is valid first !
00385         if(!(*hit).isValid()) {
00386           LogDebug("") << " InValid hit skipped in coarseHitSelection " << std::endl ;
00387           continue ;
00388         }
00389         std::string theDet = "null";
00390         int theLayer = -999;
00391         bool isStereoDet = false ;
00392         if(tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TIB) { 
00393           theDet = "TIB" ;
00394           theLayer = TIBDetId(id).layer(); 
00395           if(TIBDetId(id).stereo()==1) { isStereoDet = true ; }
00396         } else if
00397           (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TOB) { 
00398           theDet = "TOB" ;
00399           theLayer = TOBDetId(id).layer(); 
00400           if(TOBDetId(id).stereo()==1) { isStereoDet = true ; }
00401         }else if
00402           (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TID) { 
00403           theDet = "TID" ;
00404           theLayer = TIDDetId(id).wheel();  // or ring  ?
00405           if(TIDDetId(id).stereo()==1) { isStereoDet = true ; }
00406         }else if
00407           (tracker_p_->idToDetUnit(hit->geographicalId())->type().subDetector() == GeomDetEnumerators::TEC) { 
00408           theDet = "TEC" ;
00409           theLayer = TECDetId(id).wheel();  // or ring or petal ?
00410           if(TECDetId(id).stereo()==1) { isStereoDet = true ; }
00411         } else {
00412           LogDebug("") << " UHOH BIG PROBLEM - Unrecognized SI Layer" ;
00413           LogDebug("") << " Det "<< theDet << " Lay " << theLayer ;
00414           assert(1!=1) ;
00415         }
00416 
00417         if ((endcap  &&  stereo && (theDet=="TID" || theDet== "TEC") && isStereoDet ) ||
00418             (endcap  &&  !stereo && (theDet=="TID" || theDet== "TEC") && !isStereoDet )  ||
00419             (!endcap  && stereo && (theDet=="TIB" || theDet=="TOB") &&  isStereoDet )    ||
00420             (!endcap  &&  !stereo && (theDet=="TIB" || theDet=="TOB" )&& !isStereoDet )
00421             ) {  
00422               
00423 
00424           hitPointersOut.push_back(&(*hit));
00425 
00426         } // end if this is the right subdetector
00427       } // end loop over hits
00428     } // end if this detector id doesn't have too many hits on it
00429   } // end loop over detector ids
00430 }
00431 
00432 // select all matched hits for now
00433 
00434 void SiStripElectronAlgo::coarseMatchedHitSelection(std::vector<const SiStripMatchedRecHit2D*>& coarseMatchedHitPointersOut)
00435 {
00436   
00437   // Loop over the detector ids
00438   SiStripMatchedRecHit2DCollection::const_iterator itdet = matchedHits_p_->begin(), eddet = matchedHits_p_->end();
00439   for (; itdet != eddet; ++itdet) {
00440     
00441     // Get the hits on this detector id
00442     SiStripMatchedRecHit2DCollection::DetSet hits = *itdet ;
00443     
00444     // Count the number of hits on this detector id
00445     unsigned int numberOfHits = 0;
00446     for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = hits.begin();  hit != hits.end();  ++hit) {
00447       if ( !((hit->geographicalId()).subdetId() == StripSubdetector::TIB) &&
00448            !( (hit->geographicalId()).subdetId() == StripSubdetector::TOB )) { break;}
00449       numberOfHits++;
00450       if (numberOfHits > maxHitsOnDetId_) { break; }
00451     }
00452     
00453     // Only take the hits if there aren't too many
00454     if (numberOfHits <= maxHitsOnDetId_) {
00455       for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = hits.begin();  hit != hits.end();  ++hit) {
00456         if(!(*hit).isValid()) {
00457           LogDebug("") << " InValid hit skipped in coarseMatchedHitSelection " << std::endl ;
00458           continue ;
00459         }
00460         if ( !((hit->geographicalId()).subdetId() == StripSubdetector::TIB) &&
00461              !( (hit->geographicalId()).subdetId() == StripSubdetector::TOB )) { break;}
00462         
00463         coarseMatchedHitPointersOut.push_back(&(*hit));
00464       } // end loop over hits
00465       
00466     } // end if this detector id doesn't have too many hits on it
00467   } // end loop over detector ids
00468   
00469 
00470 }// end of matchedHitSelection
00471 
00472 
00473 
00474 
00475 // projects a phi band of width phiBandWidth_ from supercluster into tracker (given a chargeHypothesis)
00476 // fills *_pos_ or *_neg_ member data with the results
00477 // returns true iff the electron/positron passes cuts
00478 bool SiStripElectronAlgo::projectPhiBand(float chargeHypothesis, const reco::SuperClusterRef& superclusterIn)
00479 {
00480   // This algorithm projects a phi band into the tracker three times:
00481   // (a) for all stereo hits, (b) for barrel rphi hits, and (c) for
00482   // endcap zphi hits.  While accumulating stereo hits in step (a),
00483   // we fit r vs z to a line.  This resolves the ambiguity in z for
00484   // rphi hits and the ambiguity in r for zphi hits.  We can then cut
00485   // on the z of rphi hits (a little wider than one strip length),
00486   // and we can convert the z of zphi hits into r to apply the phi
00487   // band cut.  (We don't take advantage of the endcap strips'
00488   // segmentation in r.)
00489   // 
00490   // As we project a phi band into the tracker, we count hits within
00491   // that band and performs a linear fit for phi vs r.  The number of
00492   // hits and reduced chi^2 from the fit are used to select a good
00493   // candidate.
00494 
00495   // set limits on Si hits - smallest error that makes sense = 1 micron ?
00496   static const double rphiHitSmallestError = 0.0001 ;
00497   static const double stereoHitSmallestError = 0.0001 ;
00498   //not used since endcap is not implemented  
00499   // static const double zphiHitSmallestError = 0.0001 ;
00500 
00501 
00502   static const double stereoPitchAngle = 0.100 ; // stereo angle in rad
00503   static const double cos2SiPitchAngle = cos(stereoPitchAngle)*cos(stereoPitchAngle) ;
00504   static const double sin2SiPitchAngle = sin(stereoPitchAngle)*sin(stereoPitchAngle) ;
00505   // overall misalignment fudge to be added in quad to position errors.
00506   // this is a rough approx to values reported in tracking meet 5/16/2007
00507   static const double rphiErrFudge = 0.0200 ;
00508   static const double stereoErrFudge = 0.0200 ;
00509 
00510   // max chi2 of a hit on an SiDet relative to the prediction
00511   static const double chi2HitMax = 25.0 ;
00512 
00513   // Minimum number of hits to consider on a candidate
00514   static const unsigned int nHitsLeftMinimum = 3  ;
00515 
00516   // Create and fill vectors of pointers to hits
00517   std::vector<const SiStripRecHit2D*> stereoHits;
00518   std::vector<const SiStripRecHit2D*> rphiBarrelHits;
00519   std::vector<const SiStripRecHit2D*> zphiEndcapHits;
00520 
00521   //                                 stereo? endcap?
00522   coarseHitSelection(stereoHits,     true,   false);
00523 
00524   // skip endcap stereo for now
00525   //  LogDebug("") << " Getting endcap stereo hits " ;
00526   // coarseHitSelection(stereoHits,     true,   true);
00527 
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       OmniClusterRef const & stereocluster=(*hit)->omniClusterRef();
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             OmniClusterRef const &  mystereocluster = matchedhitlist[i]->stereoClusterRef();
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       OmniClusterRef const & monocluster=(*hit)->omniClusterRef();
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               OmniClusterRef const &  mymonocluster = matchedhitlist[i]->monoClusterRef();
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 #ifdef EDM_ML_DEBUG 
00823   std::ostringstream debugstr5;
00824   debugstr5 << " List of hits before uniqification " << " \n" ;
00825   for (unsigned int i = 0;  i < uselist.size();  i++) {
00826     if ( uselist[i] ) {
00827       const SiStripRecHit2D* hit = hitlist[i];
00828       debugstr5 << " DetID " << ((hit)->geographicalId()).rawId()
00829                 << " R " << rlist[i] 
00830                 << " Phi " << philist[i]
00831                 << " Weight " << w2list[i] 
00832                 << " PhiPred " << (rlist[i]-scr)*phiVsRSlope  
00833                 << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
00834                 << " \n" ;
00835     }
00836   }
00837   debugstr5 << " \n\n\n" ;
00838   
00839   debugstr5 << " Counting number of unique detectors " << " \n" ;
00840 
00841   debugstr5 << " These are all the detectors with hits " << " \n";
00842 #endif
00843   // Loop over hits, and find the best hit for each SiDetUnit - keep only those
00844   // get a list of DetIds in uselist
00845   std::vector<uint32_t> detIdList ;
00846 
00847   for (unsigned int i = 0;  i < uselist.size();  i++) {
00848     if (uselist[i]) {
00849       const SiStripRecHit2D* hit = hitlist[i];
00850       uint32_t detIDUsed = ((hit)->geographicalId()).rawId() ;
00851 #ifdef EDM_ML_DEBUG 
00852       debugstr5<< " DetId " << detIDUsed << "\n";
00853 #endif
00854       detIdList.push_back(detIDUsed) ;
00855     }
00856   }
00857 
00858 #ifdef EDM_ML_DEBUG 
00859   debugstr5 << " There are " << detIdList.size() << " hits on detectors \n";
00860 #endif
00861   // now sort and then uniq this list of detId
00862   std::sort(detIdList.begin(), detIdList.end());
00863   detIdList.erase(
00864                   std::unique(detIdList.begin(), detIdList.end()), detIdList.end());
00865 #ifdef EDM_ML_DEBUG 
00866   debugstr5 << " There are " << detIdList.size() << " unique detectors \n";
00867 #endif
00868   //now we have a list of unique detectors used
00869 
00870 
00871 #ifdef EDM_ML_DEBUG 
00872   debugstr5 << " Looping over detectors to choose best hit " << " \n";
00873 #endif
00874   // Loop over detectors ID and hits to create list of hits on same DetId
00875   for (unsigned int idet = 0 ; idet < detIdList.size() ; idet++ ) {
00876     for (unsigned int i = 0;  i+1 < uselist.size();  i++) {
00877       if (uselist[i]) {
00878         // Get Chi2 of this hit relative to predicted hit
00879         const SiStripRecHit2D* hit1 = hitlist[i];
00880         double r_hit1 = rlist[i];
00881         double phi_hit1 = philist[i];
00882         double w2_hit1 = w2list[i] ;
00883         double phi_pred1 = (r_hit1-scr)*phiVsRSlope ; 
00884         double chi1 = (phi_hit1-phi_pred1)*(phi_hit1-phi_pred1)*w2_hit1;
00885         if(detIdList[idet]== ((hit1)->geographicalId()).rawId()) {
00886           for (unsigned int j = i+1;  j < uselist.size();  j++) {
00887             if (uselist[j]) {
00888               const SiStripRecHit2D* hit2 = hitlist[j];
00889               if(detIdList[idet]== ((hit2)->geographicalId()).rawId()) {
00890 #ifdef EDM_ML_DEBUG 
00891                 debugstr5 << " Found 2 hits on same Si Detector " 
00892                           << ((hit2)->geographicalId()).rawId() << "\n";
00893 #endif
00894                 // Get Chi2 of this hit relative to predicted hit
00895                 double r_hit2 = rlist[j];
00896                 double phi_hit2 = philist[j];
00897                 double w2_hit2 = w2list[j] ;
00898                 double phi_pred2 = (r_hit2-scr)*phiVsRSlope ; 
00899                 double chi2 = (phi_hit2-phi_pred2)*(phi_hit2-phi_pred2)*w2_hit2;
00900 #ifdef EDM_ML_DEBUG 
00901                 debugstr5 << " Chi1 " << chi1 << " Chi2 " << chi2 <<"\n";
00902 #endif
00903                 if( chi1< chi2 ){
00904                   uselist[j] = 0;
00905                 }else{
00906                   uselist[i] = 0;
00907                 }
00908 
00909               } // end of Det check
00910             } // end of j- usehit check
00911           } // end of j hit list loop
00912           
00913         } // end of detector check
00914       } // end of uselist check
00915     } // end of i-hit list loop
00916 
00917 
00918   } // end of DetId loop
00919 
00920 
00921   
00922   // now let's through out hits with a predicted chi > chi2HitMax 
00923   for ( unsigned int i = 0;  i < uselist.size();  i++ ) { 
00924     if ( uselist[i] ) { 
00925       double localchi2 = (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i] ;
00926       if(localchi2 > chi2HitMax ) {
00927 #ifdef EDM_ML_DEBUG 
00928         const SiStripRecHit2D* hit = hitlist[i];
00929         debugstr5 << " Throwing out "
00930                   <<" DetID " << ((hit)->geographicalId()).rawId()
00931                   << " R " << rlist[i] 
00932                   << " Phi " << philist[i]
00933                   << " Weight " << w2list[i] 
00934                   << " PhiPred " << (rlist[i]-scr)*phiVsRSlope  
00935                   << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
00936                   << " \n" ;
00937 #endif
00938         uselist[i]=0 ;
00939       }
00940     }
00941   }
00942 
00943 #ifdef EDM_ML_DEBUG 
00944   debugstr5 << " List of hits after uniqification " << " \n" ;
00945   for (unsigned int i = 0;  i < uselist.size();  i++) {
00946     if ( uselist[i] ) {
00947       const SiStripRecHit2D* hit = hitlist[i];
00948       debugstr5 << " DetID " << ((hit)->geographicalId()).rawId()
00949                 << " R " << rlist[i] 
00950                 << " Phi " << philist[i]
00951                 << " Weight " << w2list[i] 
00952                 << " PhiPred " << (rlist[i]-scr)*phiVsRSlope  
00953                 << " Chi2 " << (philist[i]-(rlist[i]-scr)*phiVsRSlope)*(philist[i]-(rlist[i]-scr)*phiVsRSlope)*w2list[i]
00954                 << " \n" ;
00955     }
00956   }
00957   debugstr5 << " \n\n\n" ;
00958 #endif
00959 
00960   // need to check that there are more than 2 hits left here!
00961   unsigned int nHitsLeft =0;
00962   for (unsigned int i = 0;  i < uselist.size();  i++) {
00963     if ( uselist[i] ) {
00964       nHitsLeft++;
00965     }
00966   }
00967   if(nHitsLeft < nHitsLeftMinimum ) {
00968 #ifdef EDM_ML_DEBUG 
00969     debugstr5 << " Too few hits "<< nHitsLeft 
00970               << " left - return false " << " \n";
00971 #endif
00972     return false ;
00973   }
00974 #ifdef EDM_ML_DEBUG 
00975   LogDebug("") << debugstr5.str();
00976 #endif
00977 
00978   
00979   // Calculate a linear phi(r) fit and drop hits until the biggest contributor to chi^2 is less than maxNormResid_
00980   bool done = false;
00981   double intercept = 0 ;
00982   double slope = 0 ;
00983   double chi2 = 0;
00984 
00985   std::ostringstream debugstr4;
00986   debugstr4 <<" Calc of phi(r) "<<" \n";
00987 
00988   while (!done) {
00989     // Use an iterative update of <psi>, <r> etc instead in an
00990     // attempt to minize divisions of  large numbers
00991     // The linear fit  psi = slope*r + intercept
00992     //             slope is (<psi*r>-<psi>*<r>) / (<r^2>-<r>^2>)
00993     //             intercept is <psi> - slope*<r>
00994     
00995     double phiBar   = 0.;
00996     double phiBarOld   = 0.;
00997     double rBar  = 0.;
00998     double rBarOld  = 0.;
00999     double r2Bar = 0.;
01000     double r2BarOld = 0.;
01001     double phirBar = 0.;
01002     double phirBarOld = 0.;
01003     double totalWeight = 0.;
01004     double totalWeightOld = 0.; 
01005     unsigned int uselist_size = uselist.size();
01006     for (unsigned int i = 0;  i < uselist_size;  i++) {
01007       if (uselist[i]) {
01008         double r = rlist[i];
01009         double phi = philist[i];
01010         double weight2 = w2list[i];
01011         
01012         totalWeight = totalWeightOld + weight2 ;
01013         
01014         //weight2 is 1/sigma^2. Typically sigma is 100micron pitch
01015         // over root(12) = 30 microns so weight2 is about 10^5-10^6
01016         
01017         double totalWeightRatio = totalWeightOld/totalWeight ;
01018         double localWeightRatio = weight2/totalWeight ;
01019         
01020         phiBar = phiBarOld*totalWeightRatio + phi*localWeightRatio ; 
01021         rBar = rBarOld*totalWeightRatio + r*localWeightRatio ; 
01022         r2Bar= r2BarOld*totalWeightRatio + r*r*localWeightRatio ; 
01023         phirBar = phirBarOld*totalWeightRatio + phi*r*localWeightRatio ;
01024         
01025         totalWeightOld = totalWeight ;
01026         phiBarOld = phiBar ;
01027         rBarOld = rBar ;
01028         r2BarOld = r2Bar ;
01029         phirBarOld = phirBar ;  
01030 #ifdef EDM_ML_DEBUG 
01031         debugstr4 << " totalWeight " << totalWeight 
01032                   << " totalWeightRatio " << totalWeightRatio
01033                   << " localWeightRatio "<< localWeightRatio
01034                   << " phiBar " << phiBar
01035                   << " rBar " << rBar 
01036                   << " r2Bar " << r2Bar
01037                   << " phirBar " << phirBar 
01038                   << " \n ";
01039 #endif
01040 
01041       } // end of use hit loop
01042     } // end of hit loop to calculate a linear fit
01043     slope = (phirBar-phiBar*rBar)/(r2Bar-rBar*rBar);
01044     intercept = phiBar-slope*rBar ;
01045 
01046     debugstr4 << " end of loop slope " << slope 
01047               << " intercept " << intercept << " \n";
01048 
01049     // Calculate chi^2
01050     chi2 = 0.;
01051     double biggest_normresid = -1.;
01052     unsigned int biggest_index = 0;
01053     for (unsigned int i = 0;  i < uselist_size;  i++) {
01054       if (uselist[i]) {
01055         double r = rlist[i];
01056         double phi = philist[i];
01057         double weight2 = w2list[i];
01058 
01059         double normresid = weight2 * (phi - intercept - slope*r)*(phi - intercept - slope*r);
01060         chi2 += normresid;
01061 
01062         if (normresid > biggest_normresid) {
01063           biggest_normresid = normresid;
01064           biggest_index = i;
01065         }
01066       }
01067     } // end loop over hits to calculate chi^2 and find its biggest contributer
01068 
01069     if (biggest_normresid > maxNormResid_) {
01070 #ifdef EDM_ML_DEBUG
01071       debugstr4 << "Dropping hit from fit due to Chi2 " << " \n" ;
01072       const SiStripRecHit2D* hit = hitlist[biggest_index];
01073       debugstr4 << " DetID " << ((hit)->geographicalId()).rawId()
01074                 << " R " << rlist[biggest_index]
01075                 << " Phi " << philist[biggest_index]
01076                 << " Weight " << w2list[biggest_index]
01077                 << " PhiPred " << (rlist[biggest_index]-scr)*phiVsRSlope  
01078                 << " Chi2 " << (philist[biggest_index]-(rlist[biggest_index]-scr)*phiVsRSlope)*(philist[biggest_index]-(rlist[biggest_index]-scr)*phiVsRSlope)*w2list[biggest_index] 
01079                 << " normresid " <<  biggest_normresid
01080                 << " \n\n";
01081 #endif
01082       uselist[biggest_index] = false;
01083     }
01084     else {
01085       done = true;
01086     }
01087   } // end loop over trial fits
01088 #ifdef EDM_ML_DEBUG 
01089   debugstr4 <<" Fit " 
01090             << " Intercept  " << intercept
01091             << " slope " << slope 
01092             << " chi2 " << chi2
01093             << " \n" ;
01094 
01095   LogDebug("") <<   debugstr4.str() ;
01096 #endif
01097 
01098   // Now we have intercept, slope, and chi2; uselist to tell us which hits are used, and hitlist for the hits
01099 
01100   // Identify the innermost hit
01101   const SiStripRecHit2D* innerhit = (SiStripRecHit2D*)(0);
01102   double innerhitRadius = -1.;  // meaningless until innerhit is defined
01103 
01104   // Copy hits into an OwnVector, which we put in the TrackCandidate
01105   std::vector<const TrackingRecHit*> outputHits;
01106   // Reference rphi and stereo hits into RefVectors, which we put in the SiStripElectron
01107   std::vector<SiStripRecHit2D> outputRphiHits;
01108   std::vector<SiStripRecHit2D> outputStereoHits;
01109 
01110   typedef edm::Ref<SiStripRecHit2DCollection,SiStripRecHit2D> SiStripRecHit2DRef;
01111 
01112 
01113   for (unsigned int i = 0;  i < uselist.size();  i++) {
01114     if (uselist[i]) {
01115       double r = rlist[i];
01116       const SiStripRecHit2D* hit = hitlist[i];
01117 
01118       // Keep track of the innermost hit
01119       if (innerhit == (SiStripRecHit2D*)(0)  ||  r < innerhitRadius) {
01120         innerhit = hit;
01121         innerhitRadius = r;
01122       }
01123          
01124       if (typelist[i] == 0) {
01125         numberOfStereoHits++;
01126 
01127         // Copy this hit for the TrajectorySeed
01128         outputHits.push_back(hit);
01129         outputStereoHits.push_back(*hit);
01130       }
01131       else if (typelist[i] == 1) {
01132         numberOfBarrelRphiHits++;
01133 
01134         // Copy this hit for the TrajectorySeed
01135         outputHits.push_back(hit);
01136         outputRphiHits.push_back(*hit);
01137       }
01138       else if (typelist[i] == 2) {
01139         numberOfEndcapZphiHits++;
01140 
01141         // Copy this hit for the TrajectorySeed
01142         outputHits.push_back(hit);
01143         outputRphiHits.push_back(*hit);
01144       }
01145     }
01146   } // end loop over all hits, after having culled the ones with big residuals
01147 
01148   unsigned int totalNumberOfHits = numberOfStereoHits + numberOfBarrelRphiHits + numberOfEndcapZphiHits;
01149   double reducedChi2 = (totalNumberOfHits > 2 ? chi2 / (totalNumberOfHits - 2) : 1e10);
01150 
01151   // Select this candidate if it passes minHits_ and maxReducedChi2_ cuts
01152   if (totalNumberOfHits >= minHits_  &&  reducedChi2 <= maxReducedChi2_) {
01153     // GlobalTrajectoryParameters evaluated at the position of the innerhit
01154     GlobalPoint position = tracker_p_->idToDet(innerhit->geographicalId())->surface().toGlobal(innerhit->localPosition());
01155 
01156     // Use our phi(r) linear fit to correct pT (pT is inversely proportional to slope)
01157     // (By applying a correction instead of going back to first
01158     // principles, any correction to the phiVsRSlope definition
01159     // above will be automatically propagated here.)
01160     double correct_pT = pT * phiVsRSlope / slope;
01161 
01162     // Our phi(r) linear fit returns phi relative to the supercluster phi for a given radius
01163     double phi = intercept + slope*sqrt(position.x()*position.x() + position.y()*position.y()) + superclusterIn->position().phi();
01164 
01165     // Our z(r) linear fit gives us a better measurement of eta/dip angle
01166     double pZ = correct_pT * zVsRSlope;
01167 
01168     // Now we can construct a trajectory momentum out of linear fits to hits
01169     GlobalVector momentum = GlobalVector(correct_pT*cos(phi), correct_pT*sin(phi), pZ);
01170 
01171     if (chargeHypothesis > 0.) {
01172       redchi2_pos_ = chi2 / double(totalNumberOfHits - 2);
01173       position_pos_ = position;
01174       momentum_pos_ = momentum;
01175       innerhit_pos_ = innerhit;
01176       outputHits_pos_ = outputHits;
01177       outputRphiHits_pos_ = outputRphiHits;
01178       outputStereoHits_pos_ = outputStereoHits;
01179       phiVsRSlope_pos_ = phiVsRSlope;
01180       slope_pos_ = slope;
01181       intercept_pos_ = intercept;
01182       chi2_pos_ = chi2;
01183       ndof_pos_ = totalNumberOfHits - 2;
01184       correct_pT_pos_ = correct_pT;
01185       pZ_pos_ = pZ;
01186       zVsRSlope_pos_ = zVsRSlope;
01187       numberOfStereoHits_pos_ = numberOfStereoHits;
01188       numberOfBarrelRphiHits_pos_ = numberOfBarrelRphiHits;
01189       numberOfEndcapZphiHits_pos_ = numberOfEndcapZphiHits;
01190     }
01191     else {
01192       redchi2_neg_ = chi2 / double(totalNumberOfHits - 2);
01193       position_neg_ = position;
01194       momentum_neg_ = momentum;
01195       innerhit_neg_ = innerhit;
01196       outputHits_neg_ = outputHits;
01197       outputRphiHits_neg_ = outputRphiHits;
01198       outputStereoHits_neg_ = outputStereoHits;
01199       phiVsRSlope_neg_ = phiVsRSlope;
01200       slope_neg_ = slope;
01201       intercept_neg_ = intercept;
01202       chi2_neg_ = chi2;
01203       ndof_neg_ = totalNumberOfHits - 2;
01204       correct_pT_neg_ = correct_pT;
01205       pZ_neg_ = pZ;
01206       zVsRSlope_neg_ = zVsRSlope;
01207       numberOfStereoHits_neg_ = numberOfStereoHits;
01208       numberOfBarrelRphiHits_neg_ = numberOfBarrelRphiHits;
01209       numberOfEndcapZphiHits_neg_ = numberOfEndcapZphiHits;
01210     }
01211 
01212     LogDebug("") << " return from projectFindBand with True \n" << std::endl ;
01213     return true;
01214   } // end if this is a good electron candidate
01215 
01216   // Signal for a failed electron candidate
01217   LogDebug("") << " return from projectFindBand with False \n" << std::endl ;
01218   return false;
01219 }
01220 
01221 //
01222 // const member functions
01223 //
01224 
01225 //
01226 // static member functions
01227 //
01228