CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

TrackerSeedCleaner Class Reference

#include <TrackerSeedCleaner.h>

List of all members.

Public Types

typedef std::vector
< TrajectorySeed
tkSeeds

Public Member Functions

virtual void clean (const reco::TrackRef &, const RectangularEtaPhiTrackingRegion &region, tkSeeds &)
 the cleaner
virtual void init (const MuonServiceProxy *service)
 intizialization
virtual void setEvent (const edm::Event &)
 setEvent
 TrackerSeedCleaner (const edm::ParameterSet &pset)
 constructor
virtual ~TrackerSeedCleaner ()
 destructor

Private Attributes

edm::Handle< reco::BeamSpotbsHandle_
std::string builderName_
bool cleanBySharedHits
edm::InputTag theBeamSpotTag
const edm::EventtheEvent
const MuonServiceProxytheProxyService
RedundantSeedCleanertheRedundantCleaner
edm::ESHandle
< TransientTrackingRecHitBuilder
theTTRHBuilder
bool useDirection_Cleaner
bool usePt_Cleaner

Detailed Description

Seeds Cleaner based on direction

Date:
2010/02/16 17:08:46
Revision:
1.4
Author:
A. Grelli - Purdue University, Pavia University

Definition at line 35 of file TrackerSeedCleaner.h.


Member Typedef Documentation

Definition at line 39 of file TrackerSeedCleaner.h.


Constructor & Destructor Documentation

TrackerSeedCleaner::TrackerSeedCleaner ( const edm::ParameterSet pset) [inline]

constructor

Definition at line 41 of file TrackerSeedCleaner.h.

References builderName_, cleanBySharedHits, edm::ParameterSet::getParameter(), theBeamSpotTag, useDirection_Cleaner, and usePt_Cleaner.

                                                  : theProxyService(0),theEvent(0) {
                   builderName_ = pset.getParameter<std::string>("TTRHBuilder");
                   theBeamSpotTag = pset.getParameter<edm::InputTag>("beamSpot");
                   useDirection_Cleaner = pset.getParameter<bool>("directionCleaner");
                   usePt_Cleaner = pset.getParameter<bool>("ptCleaner");
                   cleanBySharedHits = pset.getParameter<bool>("cleanerFromSharedHits");
  }
virtual TrackerSeedCleaner::~TrackerSeedCleaner ( ) [inline, virtual]

destructor

Definition at line 53 of file TrackerSeedCleaner.h.

{}

Member Function Documentation

void TrackerSeedCleaner::clean ( const reco::TrackRef muR,
const RectangularEtaPhiTrackingRegion region,
tkSeeds seeds 
) [virtual]

the cleaner

Definition at line 63 of file TrackerSeedCleaner.cc.

References Geom::deltaPhi(), TrackingRegionBase::direction(), PV3DBase< T, PVType, FrameType >::eta(), MuonErrorMatrixValues_cff::etaRange, RectangularEtaPhiTrackingRegion::etaRange(), TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateClosestToBeamLine::isValid(), LogDebug, FreeTrajectoryState::momentum(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), RectangularEtaPhiTrackingRegion::phiMargin(), FreeTrajectoryState::position(), TrackingRegionBase::ptMin(), PtMinSelector_cfg::ptMin, query::result, TrajectoryStateClosestToBeamLine::trackStateAtPCA(), TrajectoryStateTransform::transientState(), PV3DBase< T, PVType, FrameType >::x(), reco::BeamSpot::x0(), PV3DBase< T, PVType, FrameType >::y(), reco::BeamSpot::y0(), PV3DBase< T, PVType, FrameType >::z(), and reco::BeamSpot::z0().

Referenced by TSGFromL2Muon::produce().

                                                                                                                       {


 // call the shared input cleaner
 if(cleanBySharedHits) theRedundantCleaner->define(seeds);

 theProxyService->eventSetup().get<TransientRecHitRecord>().get(builderName_,theTTRHBuilder);

 LogDebug("TrackerSeedCleaner")<<seeds.size()<<" trajectory seeds to the events before cleaning"<<endl; 

 //check the validity otherwise vertexing
 const reco::BeamSpot & bs = *bsHandle_;
 /*reco track and seeds as arguments. Seeds eta and phi are checked and 
   based on deviation from L2 eta and phi seed is accepted or not*/  

 std::vector<TrajectorySeed > result;

 TrajectoryStateTransform tsTransform;
 TSCBLBuilderNoMaterial tscblBuilder;
 // PerigeeConversions tspConverter;
 for(TrajectorySeedCollection::iterator seed = seeds.begin(); seed<seeds.end(); ++seed){
        if(seed->nHits() < 2) continue; 
        //get parameters and errors from the seed state
        TransientTrackingRecHit::RecHitPointer recHit = theTTRHBuilder->build(&*(seed->recHits().second-1));
        TrajectoryStateOnSurface state = tsTransform.transientState( seed->startingState(), recHit->surface(), theProxyService->magneticField().product());

        TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed = tscblBuilder(*state.freeState(),bs);//as in TrackProducerAlgorithms
        if (!tsAtClosestApproachSeed.isValid()) continue;
        GlobalPoint vSeed1 = tsAtClosestApproachSeed.trackStateAtPCA().position();
        GlobalVector pSeed = tsAtClosestApproachSeed.trackStateAtPCA().momentum();
        GlobalPoint vSeed(vSeed1.x()-bs.x0(),vSeed1.y()-bs.y0(),vSeed1.z()-bs.z0());


        //eta,phi info from seeds 
        double etaSeed = state.globalMomentum().eta();
        double phiSeed = pSeed.phi(); 

        //if the limits are too stringent rescale limits
        typedef PixelRecoRange< float > Range;
        typedef TkTrackingRegionsMargin< float > Margin;

        Range etaRange   = region.etaRange();
        double etaLimit  = (fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) <0.1) ? 0.1 : fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) ;

        Margin phiMargin = region.phiMargin();
        double phiLimit  = (phiMargin.right() < 0.1 ) ? 0.1 : phiMargin.right(); 

        double ptSeed  = pSeed.perp();
        double ptMin   = (region.ptMin()>3.5) ? 3.5: region.ptMin();
        // Clean  
        bool inEtaRange = etaSeed >= (etaRange.mean() - etaLimit) && etaSeed <= (etaRange.mean() + etaLimit) ;
        bool inPhiRange = (fabs(deltaPhi(phiSeed,double(region.direction().phi()))) < phiLimit );
        // pt cleaner
        bool inPtRange = ptSeed >= ptMin &&  ptSeed<= 2*(muR->pt());
        
        // save efficiency don't clean triplets with pt cleaner 
        if(seed->nHits()==3) inPtRange = true;

        // use pt and angle cleaners
        if(inPtRange  && usePt_Cleaner && !useDirection_Cleaner) {

            result.push_back(*seed);
            LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed pt selection";
        }
        
        // use only angle default option
        if( inEtaRange && inPhiRange && !usePt_Cleaner && useDirection_Cleaner) {

            result.push_back(*seed);
            LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed direction selection";

        }

        // use all the cleaners
        if( inEtaRange && inPhiRange && inPtRange && usePt_Cleaner && useDirection_Cleaner) {

            result.push_back(*seed);
            LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed direction and pt selection";

        }

        
        LogDebug("TrackerSeedCleaner")<<" eta for current seed "<<etaSeed<<"\n"
                                      <<" phi for current seed "<<phiSeed<<"\n"
                                      <<" eta for L2 track  "<<muR->eta()<<"\n"
                                      <<" phi for L2 track  "<<muR->phi()<<"\n";


  }

   //the new seeds collection
   if(result.size()!=0 && (useDirection_Cleaner || usePt_Cleaner)) seeds.swap(result);

   LogDebug("TrackerSeedCleaner")<<seeds.size()<<" trajectory seeds to the events after cleaning"<<endl;
 
   return;

}
void TrackerSeedCleaner::init ( const MuonServiceProxy *  service) [virtual]

intizialization

Definition at line 45 of file TrackerSeedCleaner.cc.

Referenced by TSGFromL2Muon::beginRun().

void TrackerSeedCleaner::setEvent ( const edm::Event event) [virtual]

setEvent

Definition at line 55 of file TrackerSeedCleaner.cc.

Referenced by TSGFromL2Muon::produce().

{
 event.getByLabel(theBeamSpotTag, bsHandle_);
}

Member Data Documentation

Definition at line 65 of file TrackerSeedCleaner.h.

std::string TrackerSeedCleaner::builderName_ [private]

Definition at line 69 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().

Definition at line 64 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().

Definition at line 62 of file TrackerSeedCleaner.h.

Definition at line 61 of file TrackerSeedCleaner.h.

Definition at line 67 of file TrackerSeedCleaner.h.

Definition at line 70 of file TrackerSeedCleaner.h.

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().