CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoParticleFlow/PFTracking/src/PFDisplacedVertexCandidateFinder.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFTracking/interface/PFDisplacedVertexCandidateFinder.h"
00002 
00003 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00004 #include "DataFormats/TrackReco/interface/Track.h"
00005 #include "DataFormats/VertexReco/interface/Vertex.h"
00006 
00007 #include "MagneticField/Engine/interface/MagneticField.h"
00008 
00009 using namespace std;
00010 using namespace reco;
00011 
00012 //for debug only 
00013 //#define PFLOW_DEBUG
00014 
00015 PFDisplacedVertexCandidateFinder::PFDisplacedVertexCandidateFinder() : 
00016   vertexCandidates_( new PFDisplacedVertexCandidateCollection ),
00017   dcaCut_(1000),
00018   primaryVertexCut2_(0.0),
00019   dcaPInnerHitCut2_(1000.0),
00020   vertexCandidatesSize_(50),
00021   debug_(false) {
00022 
00023   TwoTrackMinimumDistance theMinimum(TwoTrackMinimumDistance::SlowMode);
00024   theMinimum_ = theMinimum;
00025 }
00026 
00027 
00028 PFDisplacedVertexCandidateFinder::~PFDisplacedVertexCandidateFinder() {
00029 
00030 #ifdef PFLOW_DEBUG
00031   if(debug_)
00032     cout<<"~PFDisplacedVertexCandidateFinder - number of remaining elements: "
00033         <<eventTracks_.size()<<endl;
00034 #endif
00035   
00036 }
00037 
00038 
00039 void PFDisplacedVertexCandidateFinder::setPrimaryVertex(
00040                                                edm::Handle< reco::VertexCollection > mainVertexHandle, 
00041                                                edm::Handle< reco::BeamSpot > beamSpotHandle){
00042 
00043   const math::XYZPoint beamSpot = beamSpotHandle.isValid() ? 
00044     math::XYZPoint(beamSpotHandle->x0(), beamSpotHandle->y0(), beamSpotHandle->z0()) : 
00045     math::XYZPoint(0, 0, 0);
00046 
00047   // The primary vertex is taken from the refitted list, 
00048   // if does not exist from the average offline beam spot position  
00049   // if does not exist (0,0,0) is used
00050   pvtx_ = mainVertexHandle.isValid() ? 
00051     math::XYZPoint(mainVertexHandle->begin()->x(), 
00052                    mainVertexHandle->begin()->y(), 
00053                    mainVertexHandle->begin()->z()) :
00054     beamSpot;
00055 }
00056 
00057 // Set the imput collection of tracks and calculate their
00058 // trajectory parameters the Global Trajectory Parameters
00059 void
00060 PFDisplacedVertexCandidateFinder::setInput(const edm::Handle<TrackCollection>& trackh,
00061                                            const MagneticField* magField) {
00062 
00063   magField_ = magField;
00064   trackMask_.clear();
00065   trackMask_.resize(trackh->size());
00066   eventTrackTrajectories_.clear();
00067   eventTrackTrajectories_.resize(trackh->size());
00068 
00069   for(unsigned i=0;i<trackMask_.size(); i++) trackMask_[i] = true; 
00070 
00071   eventTracks_.clear();
00072   if(trackh.isValid()) {
00073     for(unsigned i=0;i<trackh->size(); i++) {
00074 
00075       TrackRef tref( trackh, i); 
00076       TrackBaseRef tbref(tref);
00077 
00078       if( !isSelected( tbref ) ) {
00079         trackMask_[i] = false;
00080         continue;      
00081       }
00082       
00083       const Track* trk = tref.get();
00084       eventTracks_.push_back( tbref );
00085       eventTrackTrajectories_[i] = getGlobalTrajectoryParameters(trk);
00086     }
00087   }
00088 }
00089 
00090 
00091 // -------- Main function which find vertices -------- //
00092 
00093 void 
00094 PFDisplacedVertexCandidateFinder::findDisplacedVertexCandidates() {
00095 
00096   if (debug_) cout << "========= Start Finding Displaced Vertex Candidates =========" << endl;
00097   
00098   // The vertexCandidates have not been passed to the event, and need to be cleared
00099   if(vertexCandidates_.get() )vertexCandidates_->clear();
00100   else 
00101     vertexCandidates_.reset( new PFDisplacedVertexCandidateCollection );
00102 
00103   vertexCandidates_->reserve(vertexCandidatesSize_);
00104   for(IE ie = eventTracks_.begin(); 
00105       ie != eventTracks_.end();) {
00106 
00107     // Run the recursive procedure to find all tracks link together 
00108     // In one blob called Candidate
00109 
00110     PFDisplacedVertexCandidate tempVertexCandidate;
00111 
00112     ie = associate( eventTracks_.end() , ie, tempVertexCandidate);
00113 
00114     // Build remaining links in current block
00115     if(tempVertexCandidate.isValid()) {
00116       packLinks( tempVertexCandidate );
00117       vertexCandidates_->push_back( tempVertexCandidate );
00118     }
00119   }       
00120 
00121   if(debug_) cout << "========= End Finding Displaced Vertex Candidates =========" << endl;
00122 
00123 }
00124 
00125 
00126 
00127 
00128 PFDisplacedVertexCandidateFinder::IE 
00129 PFDisplacedVertexCandidateFinder::associate(IE last, IE next, 
00130                                             PFDisplacedVertexCandidate& tempVertexCandidate) {
00131 
00132     
00133 #ifdef PFLOW_DEBUG
00134   if(debug_ ) cout<<"== Start the association procedure =="<<endl;
00135 #endif
00136 
00137   if( last!= eventTracks_.end() ) {
00138     double dist = -1;
00139     GlobalPoint P(0,0,0);
00140     PFDisplacedVertexCandidate::VertexLinkTest linktest;
00141     link( (*last), (*next), dist, P, linktest); 
00142 
00143 
00144     if(dist<-0.5) {
00145 #ifdef PFLOW_DEBUG
00146       if(debug_ ) cout<<"link failed"<<endl;
00147 #endif
00148       return ++next; // association failed
00149     }
00150     else {
00151       // add next element to the current pflowblock
00152       tempVertexCandidate.addElement( (*next));  
00153       trackMask_[(*next).key()] = false;
00154 #ifdef PFLOW_DEBUG   
00155       if(debug_ ) 
00156         cout<<"link parameters "
00157             << " *next = " << (*next).key()
00158             << " *last = " << (*last).key()
00159             << "  dist = " << dist
00160             << " P.x = " << P.x() 
00161             << " P.y = " << P.y()
00162             << " P.z = " << P.z()
00163             << endl;
00164 #endif
00165     }
00166   }
00167   else {
00168     // add next element to this eflowblock
00169 #ifdef PFLOW_DEBUG   
00170     if(debug_ ) cout<<"adding to block element "
00171                     << (*next).key()
00172                     <<endl;
00173 #endif
00174     tempVertexCandidate.addElement( (*next));  
00175     trackMask_[(*next).key()] = false;
00176 
00177 
00178   }
00179 
00180   // recursive call: associate next and other unused elements 
00181 #ifdef PFLOW_DEBUG  
00182   if(debug_ ) {
00183     for(unsigned i=0; i<trackMask_.size(); i++) 
00184       cout << " Mask[" << i << "] = " << trackMask_[i];
00185     cout << "" << endl;
00186   }
00187 #endif
00188 
00189   for(IE ie = eventTracks_.begin(); 
00190       ie != eventTracks_.end();) {
00191     
00192     if( ie == last || ie == next ) { ++ie; continue;} 
00193 
00194     // *ie already included to a block
00195     if( !trackMask_[(*ie).key()] ) { ++ie; continue;}
00196 
00197 #ifdef PFLOW_DEBUG      
00198     if(debug_ ) cout<<"calling associate "
00199                     << (*next).key()
00200                     <<" & "
00201                     << (*ie).key()
00202                     <<endl;
00203 #endif
00204     ie = associate(next, ie, tempVertexCandidate);
00205     
00206   }       
00207 
00208 #ifdef PFLOW_DEBUG
00209   if(debug_ ) {
00210     cout<<"**** removing element "<<endl;
00211   }
00212 #endif
00213 
00214   IE iteratorToNextFreeElement = eventTracks_.erase( next );
00215 
00216 #ifdef PFLOW_DEBUG
00217   if(debug_ ) cout<< "== End the association procedure ==" <<endl;
00218 #endif
00219 
00220   return iteratorToNextFreeElement;
00221 }
00222 
00223 
00224 
00225 void 
00226 PFDisplacedVertexCandidateFinder::link(const TrackBaseRef& el1, 
00227                                        const TrackBaseRef& el2, 
00228                                        double& dist, GlobalPoint& P,
00229                                        PFDisplacedVertexCandidate::VertexLinkTest& vertexLinkTest
00230                                        ){
00231  
00232   if ( fabs(el1->eta()-el2->eta()) > 1) {dist = -1; return;}
00233   if ( el1->pt()>2 && el2->pt()>2 && fabs(el1->phi()-el2->phi()) > 1) {dist = -1; return;}
00234 
00235   const GlobalTrajectoryParameters& gt1 = eventTrackTrajectories_[el1.key()];
00236   const GlobalTrajectoryParameters& gt2 = eventTrackTrajectories_[el2.key()];
00237   
00238 
00239   // Closest approach
00240   theMinimum_.calculate(gt1,gt2);
00241 
00242   // Fill the parameters
00243   dist    = theMinimum_.distance();
00244   P       = theMinimum_.crossingPoint();
00245   
00246   vertexLinkTest = PFDisplacedVertexCandidate::LINKTEST_DCA; //rechit by default 
00247 
00248 
00249   // Check if the link test fails
00250   if (dist > dcaCut_) {dist = -1; return;}
00251 
00252   // Check if the closses approach point is too close to the primary vertex/beam pipe 
00253   double rho2 = P.x()*P.x()+P.y()*P.y();
00254 
00255   if ( rho2 < primaryVertexCut2_) {dist = -1; return;}
00256 
00257   // Check if the closses approach point is out of range of the tracker 
00258 
00259   double tob_rho_limit2 = 10000;
00260   double tec_z_limit = 270;
00261 
00262   if ( rho2 > tob_rho_limit2) {dist = -1; return;}
00263   if ( fabs(P.z()) > tec_z_limit) {dist = -1; return;}
00264 
00265   return;
00266 
00267 
00268   /*
00269   // Check if the inner hit of one of the tracks is not too far away from the vertex
00270   double dx = el1->innerPosition().x()-P.x();
00271   double dy = el1->innerPosition().y()-P.y();
00272   double dz = el1->innerPosition().z()-P.z();
00273   double dist2 = dx*dx+dy*dy+dz*dz; 
00274 
00275   if (dist2 > dcaPInnerHitCut2_) {
00276     dist = -1; 
00277     if (debug_) cout << "track " << el1.key() << " dist to vertex " << sqrt(dist2) << endl; 
00278     return;
00279   }
00280 
00281   dx = el2->innerPosition().x()-P.x();
00282   dy = el2->innerPosition().y()-P.y();
00283   dz = el2->innerPosition().z()-P.z();
00284   dist2 = dx*dx+dy*dy+dz*dz; 
00285 
00286   if (dist2 > dcaPInnerHitCut2_) {
00287     dist = -1; 
00288     if (debug_) cout << "track " << el2.key() << " dist to vertex " << sqrt(dist2) << endl; 
00289     return;
00290   }
00291   */
00292   
00293 }
00294 
00295 
00296 // Build up a matrix of all the links between different tracks in 
00297 // In the Candidate
00298 void 
00299 PFDisplacedVertexCandidateFinder::packLinks( PFDisplacedVertexCandidate& vertexCandidate) {
00300   
00301   
00302   const vector < TrackBaseRef >& els = vertexCandidate.elements();
00303   
00304   //First Loop: update all link data
00305   for( unsigned i1=0; i1<els.size(); i1++ ) {
00306     for( unsigned i2=i1+1; i2<els.size(); i2++ ) {
00307       
00308       // no reflexive link
00309       if( i1==i2 ) continue;
00310       
00311       double dist = -1;
00312       GlobalPoint P(0,0,0); 
00313       PFDisplacedVertexCandidate::VertexLinkTest linktest; 
00314 
00315       link( els[i1], els[i2], dist, P, linktest);
00316 
00317 
00318 #ifdef PFLOW_DEBUG
00319       if( debug_ )
00320         cout << "Setting link between elements " << i1 << " key " << els[i1].key() 
00321              << " and " << i2 << " key " << els[i2].key() 
00322              << " of dist =" << dist << " computed from link test "
00323              << linktest << endl;
00324 #endif
00325 
00326       if(dist >-0.5) vertexCandidate.setLink( i1, i2, dist, P, linktest );
00327     }
00328   }
00329 
00330 }
00331 
00332 
00333 
00334 // --------------- TOOLS -------------- //
00335 
00336 
00337 
00338 // This tool is a copy from VZeroFinder
00339 GlobalTrajectoryParameters
00340 PFDisplacedVertexCandidateFinder::getGlobalTrajectoryParameters
00341 (const Track* track) const
00342 {
00343 
00344   const GlobalPoint position(track->vx(),
00345                              track->vy(),
00346                              track->vz());
00347  
00348   const GlobalVector momentum(track->momentum().x(),
00349                               track->momentum().y(),
00350                               track->momentum().z());
00351 
00352   GlobalTrajectoryParameters gtp(position,momentum,
00353                                  track->charge(),magField_);
00354 
00355   return gtp;
00356 }
00357 
00358 
00359 // This tool allow to pre-select the track for the displaced vertex finder
00360 bool 
00361 PFDisplacedVertexCandidateFinder::goodPtResolution( const TrackBaseRef& trackref) const {
00362 
00363   double nChi2 = trackref->normalizedChi2(); 
00364   double pt = trackref->pt();
00365   double dpt = trackref->ptError();
00366   double dxy = trackref->dxy(pvtx_);
00367 
00368   double pt_error = dpt/pt*100;
00369 
00370   if (debug_) cout << " PFDisplacedVertexFinder: PFrecTrack->Track Pt= "
00371                    << pt << " dPt/Pt = " << pt_error << "% nChi2 = " << nChi2 << endl;
00372   if (nChi2 > nChi2_max_ || pt < pt_min_){
00373 
00374     if (debug_) cout << " PFBlockAlgo: skip badly measured or low pt track"
00375                      << " nChi2_cut = " << 5 
00376                      << " pt_cut = " << 0.2 << endl;
00377     return false;
00378   }
00379   //  cout << "dxy = " << dxy << endl; 
00380   if (fabs(dxy) < dxy_ && pt < pt_min_prim_) return false;
00381   //  if (fabs(dxy) < 0.2 && pt < 0.8) return false;
00382 
00383   
00384 
00385   return true;
00386 }
00387 
00388 
00389 
00390 
00391 
00392 
00393 
00394 
00395 
00396 
00397 
00398 ostream& operator<<(std::ostream& out, const PFDisplacedVertexCandidateFinder& a) {
00399   if(! out) return out;
00400   
00401   out<<"====== Particle Flow Block Algorithm ======= ";
00402   out<<endl;
00403   out<<"number of unassociated elements : "<<a.eventTracks_.size()<<endl;
00404   out<<endl;
00405   
00406   out << " Tracks selection based on " << std::endl;
00407   out << " pvtx_ = " << a.pvtx_ << std::endl;
00408   out << " fabs(dxy) < " << a.dxy_ << " and pt < "<< a.pt_min_prim_ << std::endl;
00409   out << " nChi2 < " << a.nChi2_max_ << " and pt < "<< a.pt_min_ << std::endl;
00410 
00411   out<<endl;
00412 
00413 
00414   for(PFDisplacedVertexCandidateFinder::IEC ie = a.eventTracks_.begin(); 
00415       ie != a.eventTracks_.end(); ie++) {
00416 
00417     double pt = (*ie).get()->pt(); 
00418 
00419     math::XYZPoint Pi = (*ie).get()->innerPosition(); 
00420     math::XYZPoint Po = (*ie).get()->outerPosition(); 
00421 
00422     double innermost_radius = sqrt(Pi.x()*Pi.x() + Pi.y()*Pi.y() + Pi.z()*Pi.z());
00423     double outermost_radius = sqrt(Po.x()*Po.x() + Po.y()*Po.y() + Po.z()*Po.z());
00424     double innermost_rho = sqrt(Pi.x()*Pi.x() + Pi.y()*Pi.y());
00425     double outermost_rho = sqrt(Po.x()*Po.x() + Po.y()*Po.y());
00426     
00427     out<<"ie = " << (*ie).key() 
00428        <<" pt = " << pt
00429        <<" innermost hit radius = " << innermost_radius << " rho = " << innermost_rho
00430        <<" outermost hit radius = " << outermost_radius << " rho = " << outermost_rho
00431        <<endl;
00432   }
00433 
00434 
00435   const std::auto_ptr< reco::PFDisplacedVertexCandidateCollection >& vertexCandidates
00436     = a.vertexCandidates(); 
00437     
00438   if(!vertexCandidates.get() ) {
00439     out<<"vertexCandidates already transfered"<<endl;
00440   }
00441   else {
00442     out<<"number of vertexCandidates : "<<vertexCandidates->size()<<endl;
00443     out<<endl;
00444  
00445     
00446     for(PFDisplacedVertexCandidateFinder::IBC ib=vertexCandidates->begin(); 
00447         ib != vertexCandidates->end(); ib++)
00448       ib->Dump();
00449 
00450 
00451     
00452   }
00453  
00454   return out;
00455 }