CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/RecoHI/HiTracking/plugins/HIProtoTrackSelector.h

Go to the documentation of this file.
00001 #ifndef HIProtoTrackSelection_h
00002 #define HIProtoTrackSelection_h
00003 
00004 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00005 #include "DataFormats/TrackReco/interface/Track.h"
00006 
00007 #include "DataFormats/VertexReco/interface/Vertex.h"
00008 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00009 
00010 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00011 #include "FWCore/Utilities/interface/InputTag.h"
00012 
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 
00015 #include <algorithm>
00016 #include <iostream>
00017 
00022 class HIProtoTrackSelector
00023 {
00024   
00025  public:
00026   // input collection type
00027   typedef reco::TrackCollection collection;
00028   
00029   // output collection type
00030   typedef std::vector<const reco::Track *> container;
00031   
00032   // iterator over result collection type.
00033   typedef container::const_iterator const_iterator;
00034   
00035   // constructor from parameter set configurability
00036   HIProtoTrackSelector(const edm::ParameterSet & iConfig) : 
00037     vertexCollection_(iConfig.getParameter<edm::InputTag>("VertexCollection")),
00038     beamSpotLabel_(iConfig.getParameter<edm::InputTag>("beamSpotLabel")),
00039     ptMin_(iConfig.getParameter<double>("ptMin")),      
00040     nSigmaZ_(iConfig.getParameter<double>("nSigmaZ")),
00041     minZCut_(iConfig.getParameter<double>("minZCut")),
00042     maxD0Significance_(iConfig.getParameter<double>("maxD0Significance"))
00043     {};
00044   
00045   // select object from a collection and possibly event content
00046   void select( edm::Handle<reco::TrackCollection>& TCH, const edm::Event & iEvent, const edm::EventSetup & iSetup)
00047     {
00048       selected_.clear();
00049       
00050       const collection & c = *(TCH.product());
00051       
00052       // Get fast reco vertex 
00053       edm::Handle<reco::VertexCollection> vc;
00054       iEvent.getByLabel(vertexCollection_, vc);
00055       const reco::VertexCollection * vertices = vc.product();
00056       
00057       math::XYZPoint vtxPoint(0.0,0.0,0.0);
00058       double vzErr =0.0;
00059       
00060       if(vertices->size()>0) {
00061         vtxPoint=vertices->begin()->position();
00062         vzErr=vertices->begin()->zError();
00063         edm::LogInfo("HeavyIonVertexing") << "Select prototracks compatible with median vertex"
00064                                      << "\n   vz = " << vtxPoint.Z()  
00065                                      << "\n   " << nSigmaZ_ << " vz sigmas = " << vzErr*nSigmaZ_
00066                                      << "\n   cut at = " << std::max(vzErr*nSigmaZ_,minZCut_);
00067       } 
00068       // Supress this warning, since events w/ no vertex are expected 
00069       //else {
00070         //edm::LogError("HeavyIonVertexing") << "No vertex found in collection '" << vertexCollection_ << "'";
00071       //}
00072       
00073       // Get beamspot
00074       reco::BeamSpot beamSpot;
00075       edm::Handle<reco::BeamSpot> beamSpotHandle;
00076       iEvent.getByLabel(beamSpotLabel_, beamSpotHandle);
00077       
00078       math::XYZPoint bsPoint(0.0,0.0,0.0);
00079       double bsWidth = 0.0;
00080       
00081       if ( beamSpotHandle.isValid() ) {
00082         beamSpot = *beamSpotHandle;
00083         bsPoint = beamSpot.position();
00084         bsWidth = sqrt(beamSpot.BeamWidthX()*beamSpot.BeamWidthY());
00085         edm::LogInfo("HeavyIonVertexing") << "Select prototracks compatible with beamspot"
00086                                      << "\n   (x,y,z) = (" << bsPoint.X() << "," << bsPoint.Y() << "," << bsPoint.Z() << ")"  
00087                                      << "\n   width = " << bsWidth
00088                                      << "\n   cut at d0/d0sigma = " << maxD0Significance_;
00089       } else {
00090         edm::LogError("HeavyIonVertexing") << "No beam spot available from '" << beamSpotLabel_ << "\n";
00091       }
00092       
00093       
00094       // Do selection
00095       int nSelected=0;
00096       int nRejected=0;
00097       double d0=0.0; 
00098       double d0sigma=0.0;
00099       for (reco::TrackCollection::const_iterator trk = c.begin(); trk != c.end(); ++ trk)
00100         {
00101           
00102           d0 = -1.*trk->dxy(bsPoint);
00103           d0sigma = sqrt(trk->d0Error()*trk->d0Error() + bsWidth*bsWidth);
00104           if ( trk->pt() > ptMin_ // keep only tracks above ptMin
00105                && fabs(d0/d0sigma) < maxD0Significance_ // keep only tracks with D0 significance less than cut
00106                && fabs(trk->dz(vtxPoint)) <  std::max(vzErr*nSigmaZ_,minZCut_) // within minZCut, nSigmaZ of fast vertex
00107                ) 
00108           {
00109             nSelected++;
00110             selected_.push_back( & * trk );
00111           } 
00112           else 
00113             nRejected++;
00114         }
00115       
00116       edm::LogInfo("HeavyIonVertexing") << "selected " << nSelected << " prototracks out of " << nRejected+nSelected << "\n";
00117       
00118     }
00119   
00120   // iterators over selected objects: collection begin
00121   const_iterator begin() const { return selected_.begin(); }
00122   
00123   // iterators over selected objects: collection end
00124   const_iterator end() const { return selected_.end(); }
00125   
00126   // true if no object has been selected
00127   size_t size() const { return selected_.size(); }
00128   
00129   
00130  private:
00131   container selected_;          
00132   edm::InputTag vertexCollection_; 
00133   edm::InputTag beamSpotLabel_;
00134   double ptMin_; 
00135   double nSigmaZ_;
00136   double minZCut_; 
00137   double maxD0Significance_;
00138   
00139 };
00140 
00141 
00142 #endif