CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

HIProtoTrackSelector Class Reference

#include <HIProtoTrackSelector.h>

List of all members.

Public Types

typedef reco::TrackCollection collection
typedef container::const_iterator const_iterator
typedef std::vector< const
reco::Track * > 
container

Public Member Functions

const_iterator begin () const
const_iterator end () const
 HIProtoTrackSelector (const edm::ParameterSet &iConfig)
void select (edm::Handle< reco::TrackCollection > &TCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
size_t size () const

Private Attributes

edm::InputTag beamSpotLabel_
double maxD0Significance_
double minZCut_
double nSigmaZ_
double ptMin_
container selected_
edm::InputTag vertexCollection_

Detailed Description

Selector to select prototracks that pass certain kinematic cuts based on fast vertex

Definition at line 22 of file HIProtoTrackSelector.h.


Member Typedef Documentation

Definition at line 27 of file HIProtoTrackSelector.h.

typedef container::const_iterator HIProtoTrackSelector::const_iterator

Definition at line 33 of file HIProtoTrackSelector.h.

typedef std::vector<const reco::Track *> HIProtoTrackSelector::container

Definition at line 30 of file HIProtoTrackSelector.h.


Constructor & Destructor Documentation

HIProtoTrackSelector::HIProtoTrackSelector ( const edm::ParameterSet iConfig) [inline]

Definition at line 36 of file HIProtoTrackSelector.h.

                                                        : 
    vertexCollection_(iConfig.getParameter<edm::InputTag>("VertexCollection")),
    beamSpotLabel_(iConfig.getParameter<edm::InputTag>("beamSpotLabel")),
    ptMin_(iConfig.getParameter<double>("ptMin")),      
    nSigmaZ_(iConfig.getParameter<double>("nSigmaZ")),
    minZCut_(iConfig.getParameter<double>("minZCut")),
    maxD0Significance_(iConfig.getParameter<double>("maxD0Significance"))
    {};

Member Function Documentation

const_iterator HIProtoTrackSelector::begin ( void  ) const [inline]

Definition at line 119 of file HIProtoTrackSelector.h.

References selected_.

{ return selected_.begin(); }
const_iterator HIProtoTrackSelector::end ( void  ) const [inline]

Definition at line 122 of file HIProtoTrackSelector.h.

References selected_.

{ return selected_.end(); }
void HIProtoTrackSelector::select ( edm::Handle< reco::TrackCollection > &  TCH,
const edm::Event iEvent,
const edm::EventSetup iSetup 
) [inline]

Definition at line 46 of file HIProtoTrackSelector.h.

References beamSpotLabel_, reco::BeamSpot::BeamWidthX(), reco::BeamSpot::BeamWidthY(), trackerHits::c, debug_cff::d0, edm::Event::getByLabel(), edm::HandleBase::isValid(), max(), maxD0Significance_, minZCut_, nSigmaZ_, reco::BeamSpot::position(), edm::Handle< T >::product(), ptMin_, selected_, mathSSE::sqrt(), and vertexCollection_.

    {
      selected_.clear();
      
      const collection & c = *(TCH.product());
      
      // Get fast reco vertex 
      edm::Handle<reco::VertexCollection> vc;
      iEvent.getByLabel(vertexCollection_, vc);
      const reco::VertexCollection * vertices = vc.product();
      
      math::XYZPoint vtxPoint(0.0,0.0,0.0);
      double vzErr =0.0;
      
      if(vertices->size()>0) {
        vtxPoint=vertices->begin()->position();
        vzErr=vertices->begin()->zError();
        edm::LogInfo("HeavyIonVertexing") << "Select prototracks compatible with median vertex"
                                     << "\n   vz = " << vtxPoint.Z()  
                                     << "\n   " << nSigmaZ_ << " vz sigmas = " << vzErr*nSigmaZ_
                                     << "\n   cut at = " << std::max(vzErr*nSigmaZ_,minZCut_);
      } else {
        edm::LogError("HeavyIonVertexing") << "No vertex found in collection '" << vertexCollection_ << "'";
      }
      
      // Get beamspot
      reco::BeamSpot beamSpot;
      edm::Handle<reco::BeamSpot> beamSpotHandle;
      iEvent.getByLabel(beamSpotLabel_, beamSpotHandle);
      
      math::XYZPoint bsPoint(0.0,0.0,0.0);
      double bsWidth = 0.0;
      
      if ( beamSpotHandle.isValid() ) {
        beamSpot = *beamSpotHandle;
        bsPoint = beamSpot.position();
        bsWidth = sqrt(beamSpot.BeamWidthX()*beamSpot.BeamWidthY());
        edm::LogInfo("HeavyIonVertexing") << "Select prototracks compatible with beamspot"
                                     << "\n   (x,y,z) = (" << bsPoint.X() << "," << bsPoint.Y() << "," << bsPoint.Z() << ")"  
                                     << "\n   width = " << bsWidth
                                     << "\n   cut at d0/d0sigma = " << maxD0Significance_;
      } else {
        edm::LogError("HeavyIonVertexing") << "No beam spot available from '" << beamSpotLabel_ << "\n";
      }
      
      
      // Do selection
      int nSelected=0;
      int nRejected=0;
      double d0=0.0; 
      double d0sigma=0.0;
      for (reco::TrackCollection::const_iterator trk = c.begin(); trk != c.end(); ++ trk)
        {
          
          d0 = -1.*trk->dxy(bsPoint);
          d0sigma = sqrt(trk->d0Error()*trk->d0Error() + bsWidth*bsWidth);
          if ( trk->pt() > ptMin_ // keep only tracks above ptMin
               && fabs(d0/d0sigma) < maxD0Significance_ // keep only tracks with D0 significance less than cut
               && fabs(trk->dz(vtxPoint)) <  std::max(vzErr*nSigmaZ_,minZCut_) // within minZCut, nSigmaZ of fast vertex
               ) 
          {
            nSelected++;
            selected_.push_back( & * trk );
          } 
          else 
            nRejected++;
        }
      
      edm::LogInfo("HeavyIonVertexing") << "selected " << nSelected << " prototracks out of " << nRejected+nSelected << "\n";
      
    }
size_t HIProtoTrackSelector::size ( void  ) const [inline]

Definition at line 125 of file HIProtoTrackSelector.h.

References selected_.

{ return selected_.size(); }

Member Data Documentation

Definition at line 131 of file HIProtoTrackSelector.h.

Referenced by select().

Definition at line 135 of file HIProtoTrackSelector.h.

Referenced by select().

Definition at line 134 of file HIProtoTrackSelector.h.

Referenced by select().

Definition at line 133 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::ptMin_ [private]

Definition at line 132 of file HIProtoTrackSelector.h.

Referenced by select().

Definition at line 129 of file HIProtoTrackSelector.h.

Referenced by begin(), end(), select(), and size().

Definition at line 130 of file HIProtoTrackSelector.h.

Referenced by select().