CMS 3D CMS Logo

TrackerSeedCleaner Class Reference

Seeds Cleaner based on direction
Date
2008/04/09 05:14:39
Revision
1.3
. More...

#include <RecoMuon/TrackerSeedGenerator/interface/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 MuonServiceProxy * theProxyService
RedundantSeedCleanertheRedundantCleaner
edm::ESHandle
< TransientTrackingRecHitBuilder
theTTRHBuilder
bool useDirection_Cleaner
bool usePt_Cleaner


Detailed Description

Seeds Cleaner based on direction
Date
2008/04/09 05:14:39
Revision
1.3
.

Author:
A. Grelli - Purdue University, Pavia University

Definition at line 35 of file TrackerSeedCleaner.h.


Member Typedef Documentation

typedef std::vector<TrajectorySeed> TrackerSeedCleaner::tkSeeds

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.

00041                                                   : theProxyService(0),theEvent(0) {
00042                    builderName_ = pset.getParameter<std::string>("TTRHBuilder");
00043                    theBeamSpotTag = pset.getParameter<edm::InputTag>("beamSpot");
00044                    useDirection_Cleaner = pset.getParameter<bool>("directionCleaner");
00045                    usePt_Cleaner = pset.getParameter<bool>("ptCleaner");
00046                    cleanBySharedHits = pset.getParameter<bool>("cleanerFromSharedHits");
00047   }

virtual TrackerSeedCleaner::~TrackerSeedCleaner (  )  [inline, virtual]

destructor

Definition at line 53 of file TrackerSeedCleaner.h.

00053 {}


Member Function Documentation

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

the cleaner

Definition at line 64 of file TrackerSeedCleaner.cc.

References bsHandle_, builderName_, cleanBySharedHits, RedundantSeedCleaner::define(), deltaPhi(), TrackingRegionBase::direction(), lat::endl(), PV3DBase< T, PVType, FrameType >::eta(), RectangularEtaPhiTrackingRegion::etaRange(), TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalMomentum(), LogDebug, FreeTrajectoryState::momentum(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), RectangularEtaPhiTrackingRegion::phiMargin(), FreeTrajectoryState::position(), TrackingRegionBase::ptMin(), RegionalCKFTracksForL3Isolation_cfi::ptMin, HLT_VtxMuL3::result, state, theProxyService, theRedundantCleaner, theTTRHBuilder, TrajectoryStateClosestToBeamLine::trackStateAtPCA(), TrajectoryStateTransform::transientState(), useDirection_Cleaner, usePt_Cleaner, reco::BeamSpot::x0(), reco::BeamSpot::y0(), and reco::BeamSpot::z0().

Referenced by TSGFromL2Muon::produce().

00064                                                                                                                        {
00065 
00066 
00067  // call the shared input cleaner
00068  if(cleanBySharedHits) theRedundantCleaner->define(seeds);
00069 
00070  theProxyService->eventSetup().get<TransientRecHitRecord>().get(builderName_,theTTRHBuilder);
00071 
00072  LogDebug("TrackerSeedCleaner")<<seeds.size()<<" trajectory seeds to the events before cleaning"<<endl; 
00073 
00074  //check the validity otherwise vertexing
00075  const reco::BeamSpot & bs = *bsHandle_;
00076  /*reco track and seeds as arguments. Seeds eta and phi are checked and 
00077    based on deviation from L2 eta and phi seed is accepted or not*/  
00078 
00079  std::vector<TrajectorySeed > result;
00080 
00081  TrajectoryStateTransform tsTransform;
00082  TrajectoryStateClosestToBeamLineBuilder tscblBuilder;
00083  // PerigeeConversions tspConverter;
00084  for(TrajectorySeedCollection::iterator seed = seeds.begin(); seed<seeds.end(); ++seed){
00085         if(seed->nHits() < 2) continue; 
00086         //get parameters and errors from the seed state
00087         TransientTrackingRecHit::RecHitPointer recHit = theTTRHBuilder->build(&*(seed->recHits().second-1));
00088         TrajectoryStateOnSurface state = tsTransform.transientState( seed->startingState(), recHit->surface(), theProxyService->magneticField().product());
00089 
00090         TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed = tscblBuilder(*state.freeState(),bs);//as in TrackProducerAlgorithms
00091 
00092         GlobalPoint vSeed1 = tsAtClosestApproachSeed.trackStateAtPCA().position();
00093         GlobalVector pSeed = tsAtClosestApproachSeed.trackStateAtPCA().momentum();
00094         GlobalPoint vSeed(vSeed1.x()-bs.x0(),vSeed1.y()-bs.y0(),vSeed1.z()-bs.z0());
00095 
00096 
00097         //eta,phi info from seeds 
00098         double etaSeed = state.globalMomentum().eta();
00099         double phiSeed = pSeed.phi(); 
00100 
00101         //if the limits are too stringent rescale limits
00102         typedef PixelRecoRange< float > Range;
00103         typedef TkTrackingRegionsMargin< float > Margin;
00104 
00105         Range etaRange   = region.etaRange();
00106         double etaLimit  = (fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) <0.1) ? 0.1 : fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) ;
00107 
00108         Margin phiMargin = region.phiMargin();
00109         double phiLimit  = (phiMargin.right() < 0.1 ) ? 0.1 : phiMargin.right(); 
00110 
00111         double ptSeed  = pSeed.perp();
00112         double ptMin   = (region.ptMin()>3.5) ? 3.5: region.ptMin();
00113         // Clean  
00114         bool inEtaRange = etaSeed >= (etaRange.mean() - etaLimit) && etaSeed <= (etaRange.mean() + etaLimit) ;
00115         bool inPhiRange = (fabs(deltaPhi(phiSeed,double(region.direction().phi()))) < phiLimit );
00116         // pt cleaner
00117         bool inPtRange = ptSeed >= ptMin &&  ptSeed<= 2*(muR->pt());
00118         
00119         // save efficiency don't clean triplets with pt cleaner 
00120         if(seed->nHits()==3) inPtRange = true;
00121 
00122         // use pt and angle cleaners
00123         if(inPtRange  && usePt_Cleaner && !useDirection_Cleaner) {
00124 
00125             result.push_back(*seed);
00126             LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed pt selection";
00127         }
00128         
00129         // use only angle default option
00130         if( inEtaRange && inPhiRange && !usePt_Cleaner && useDirection_Cleaner) {
00131 
00132             result.push_back(*seed);
00133             LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed direction selection";
00134 
00135         }
00136 
00137         // use all the cleaners
00138         if( inEtaRange && inPhiRange && inPtRange && usePt_Cleaner && useDirection_Cleaner) {
00139 
00140             result.push_back(*seed);
00141             LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed direction and pt selection";
00142 
00143         }
00144 
00145         
00146         LogDebug("TrackerSeedCleaner")<<" eta for current seed "<<etaSeed<<"\n"
00147                                       <<" phi for current seed "<<phiSeed<<"\n"
00148                                       <<" eta for L2 track  "<<muR->eta()<<"\n"
00149                                       <<" phi for L2 track  "<<muR->phi()<<"\n";
00150 
00151 
00152   }
00153 
00154    //the new seeds collection
00155    if(result.size()!=0 && (useDirection_Cleaner || usePt_Cleaner)) seeds.swap(result);
00156 
00157    LogDebug("TrackerSeedCleaner")<<seeds.size()<<" trajectory seeds to the events after cleaning"<<endl;
00158  
00159    return;
00160 
00161 }

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

intizialization

Definition at line 46 of file TrackerSeedCleaner.cc.

References HLT_VtxMuL3::RedundantSeedCleaner, theProxyService, and theRedundantCleaner.

Referenced by TSGFromL2Muon::beginJob().

00046                                                             {
00047 
00048   theProxyService = service;
00049   
00050   theRedundantCleaner = new RedundantSeedCleaner();
00051 }

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

setEvent

Definition at line 56 of file TrackerSeedCleaner.cc.

References bsHandle_, and theBeamSpotTag.

Referenced by TSGFromL2Muon::produce().

00057 {
00058  event.getByLabel(theBeamSpotTag, bsHandle_);
00059 }


Member Data Documentation

edm::Handle<reco::BeamSpot> TrackerSeedCleaner::bsHandle_ [private]

Definition at line 65 of file TrackerSeedCleaner.h.

Referenced by clean(), and setEvent().

std::string TrackerSeedCleaner::builderName_ [private]

Definition at line 69 of file TrackerSeedCleaner.h.

Referenced by clean(), and TrackerSeedCleaner().

bool TrackerSeedCleaner::cleanBySharedHits [private]

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by clean(), and TrackerSeedCleaner().

edm::InputTag TrackerSeedCleaner::theBeamSpotTag [private]

Definition at line 64 of file TrackerSeedCleaner.h.

Referenced by setEvent(), and TrackerSeedCleaner().

const edm::Event* TrackerSeedCleaner::theEvent [private]

Definition at line 62 of file TrackerSeedCleaner.h.

const MuonServiceProxy* TrackerSeedCleaner::theProxyService [private]

Definition at line 61 of file TrackerSeedCleaner.h.

Referenced by clean(), and init().

RedundantSeedCleaner* TrackerSeedCleaner::theRedundantCleaner [private]

Definition at line 67 of file TrackerSeedCleaner.h.

Referenced by clean(), and init().

edm::ESHandle<TransientTrackingRecHitBuilder> TrackerSeedCleaner::theTTRHBuilder [private]

Definition at line 70 of file TrackerSeedCleaner.h.

Referenced by clean().

bool TrackerSeedCleaner::useDirection_Cleaner [private]

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by clean(), and TrackerSeedCleaner().

bool TrackerSeedCleaner::usePt_Cleaner [private]

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by clean(), and TrackerSeedCleaner().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:34:10 2009 for CMSSW by  doxygen 1.5.4