CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Attributes
TrackerSeedCleaner Class Reference

#include <TrackerSeedCleaner.h>

Public Types

typedef std::vector
< TrajectorySeed
tkSeeds
 

Public Member Functions

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

Private Attributes

edm::EDGetTokenT< reco::BeamSpotbeamspotToken_
 
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

Author
A. Grelli - Purdue University, Pavia University

Definition at line 34 of file TrackerSeedCleaner.h.

Member Typedef Documentation

Definition at line 38 of file TrackerSeedCleaner.h.

Constructor & Destructor Documentation

TrackerSeedCleaner::TrackerSeedCleaner ( const edm::ParameterSet pset,
edm::ConsumesCollector iC 
)
inline

constructor

Definition at line 40 of file TrackerSeedCleaner.h.

References beamspotToken_, builderName_, cleanBySharedHits, edm::ConsumesCollector::consumes(), edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, theBeamSpotTag, useDirection_Cleaner, and usePt_Cleaner.

40  : theProxyService(0),theEvent(0) {
41  builderName_ = pset.getParameter<std::string>("TTRHBuilder");
42  theBeamSpotTag = pset.getParameter<edm::InputTag>("beamSpot");
43  useDirection_Cleaner = pset.getParameter<bool>("directionCleaner");
44  usePt_Cleaner = pset.getParameter<bool>("ptCleaner");
45  cleanBySharedHits = pset.getParameter<bool>("cleanerFromSharedHits");
47  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::InputTag theBeamSpotTag
const edm::Event * theEvent
const MuonServiceProxy * theProxyService
edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
virtual TrackerSeedCleaner::~TrackerSeedCleaner ( )
inlinevirtual

destructor

Definition at line 53 of file TrackerSeedCleaner.h.

53 {}

Member Function Documentation

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

the cleaner

Definition at line 61 of file TrackerSeedCleaner.cc.

References SiPixelRawToDigiRegional_cfi::deltaPhi, TrackingRegion::direction(), PV3DBase< T, PVType, FrameType >::eta(), MuonErrorMatrixValues_cff::etaRange, RectangularEtaPhiTrackingRegion::etaRange(), TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalMomentum(), LogDebug, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), RectangularEtaPhiTrackingRegion::phiMargin(), PtMinSelector_cfg::ptMin, TrackingRegion::ptMin(), query::result, TkTrackingRegionsMargin< T >::right(), fileCollector::seed, 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().

61  {
62 
63 
64  // call the shared input cleaner
66 
68 
69  LogDebug("TrackerSeedCleaner")<<seeds.size()<<" trajectory seeds to the events before cleaning"<<endl;
70 
71  //check the validity otherwise vertexing
72  const reco::BeamSpot & bs = *bsHandle_;
73  /*reco track and seeds as arguments. Seeds eta and phi are checked and
74  based on deviation from L2 eta and phi seed is accepted or not*/
75 
76  std::vector<TrajectorySeed > result;
77 
78 
79  TSCBLBuilderNoMaterial tscblBuilder;
80  // PerigeeConversions tspConverter;
81  for(TrajectorySeedCollection::iterator seed = seeds.begin(); seed<seeds.end(); ++seed){
82  if(seed->nHits() < 2) continue;
83  //get parameters and errors from the seed state
84  TransientTrackingRecHit::RecHitPointer recHit = theTTRHBuilder->build(&*(seed->recHits().second-1));
85  TrajectoryStateOnSurface state = trajectoryStateTransform::transientState( seed->startingState(), recHit->surface(), theProxyService->magneticField().product());
86 
87  TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed = tscblBuilder(*state.freeState(),bs);//as in TrackProducerAlgorithms
88  if (!tsAtClosestApproachSeed.isValid()) continue;
89  GlobalPoint vSeed1 = tsAtClosestApproachSeed.trackStateAtPCA().position();
90  GlobalVector pSeed = tsAtClosestApproachSeed.trackStateAtPCA().momentum();
91  GlobalPoint vSeed(vSeed1.x()-bs.x0(),vSeed1.y()-bs.y0(),vSeed1.z()-bs.z0());
92 
93 
94  //eta,phi info from seeds
95  double etaSeed = state.globalMomentum().eta();
96  double phiSeed = pSeed.phi();
97 
98  //if the limits are too stringent rescale limits
100  typedef TkTrackingRegionsMargin< float > Margin;
101 
102  Range etaRange = region.etaRange();
103  double etaLimit = (fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) <0.1) ? 0.1 : fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) ;
104 
105  Margin phiMargin = region.phiMargin();
106  double phiLimit = (phiMargin.right() < 0.1 ) ? 0.1 : phiMargin.right();
107 
108  double ptSeed = pSeed.perp();
109  double ptMin = (region.ptMin()>3.5) ? 3.5: region.ptMin();
110  // Clean
111  bool inEtaRange = etaSeed >= (etaRange.mean() - etaLimit) && etaSeed <= (etaRange.mean() + etaLimit) ;
112  bool inPhiRange = (fabs(deltaPhi(phiSeed,double(region.direction().phi()))) < phiLimit );
113  // pt cleaner
114  bool inPtRange = ptSeed >= ptMin && ptSeed<= 2*(muR->pt());
115 
116  // save efficiency don't clean triplets with pt cleaner
117  if(seed->nHits()==3) inPtRange = true;
118 
119  // use pt and angle cleaners
120  if(inPtRange && usePt_Cleaner && !useDirection_Cleaner) {
121 
122  result.push_back(*seed);
123  LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed pt selection";
124  }
125 
126  // use only angle default option
127  if( inEtaRange && inPhiRange && !usePt_Cleaner && useDirection_Cleaner) {
128 
129  result.push_back(*seed);
130  LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed direction selection";
131 
132  }
133 
134  // use all the cleaners
135  if( inEtaRange && inPhiRange && inPtRange && usePt_Cleaner && useDirection_Cleaner) {
136 
137  result.push_back(*seed);
138  LogDebug("TrackerSeedCleaner")<<" Keeping the seed : this seed passed direction and pt selection";
139 
140  }
141 
142 
143  LogDebug("TrackerSeedCleaner")<<" eta for current seed "<<etaSeed<<"\n"
144  <<" phi for current seed "<<phiSeed<<"\n"
145  <<" eta for L2 track "<<muR->eta()<<"\n"
146  <<" phi for L2 track "<<muR->phi()<<"\n";
147 
148 
149  }
150 
151  //the new seeds collection
152  if(result.size()!=0 && (useDirection_Cleaner || usePt_Cleaner)) seeds.swap(result);
153 
154  LogDebug("TrackerSeedCleaner")<<seeds.size()<<" trajectory seeds to the events after cleaning"<<endl;
155 
156  return;
157 
158 }
#define LogDebug(id)
double z0() const
z coordinate
Definition: BeamSpot.h:68
T perp() const
Definition: PV3DBase.h:72
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
void define(std::vector< TrajectorySeed > &)
collection definition
RedundantSeedCleaner * theRedundantCleaner
GlobalVector const & direction() const
the direction around which region is constructed
const MuonServiceProxy * theProxyService
FreeTrajectoryState const * freeState(bool withErrors=true) const
edm::Handle< reco::BeamSpot > bsHandle_
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
std::shared_ptr< TrackingRecHit const > RecHitPointer
PixelRecoRange< float > Range
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
edm::ESHandle< TransientTrackingRecHitBuilder > theTTRHBuilder
float ptMin() const
minimal pt of interest
T eta() const
Definition: PV3DBase.h:76
GlobalVector globalMomentum() const
double y0() const
y coordinate
Definition: BeamSpot.h:66
T x() const
Definition: PV3DBase.h:62
double x0() const
x coordinate
Definition: BeamSpot.h:64
const Range & etaRange() const
allowed eta range [eta_min, eta_max] interval
void TrackerSeedCleaner::init ( const MuonServiceProxy service)
virtual

intizialization

Definition at line 43 of file TrackerSeedCleaner.cc.

References fff_monitoring::service.

Referenced by TSGFromL2Muon::beginRun().

43  {
44 
46 
48 }
RedundantSeedCleaner * theRedundantCleaner
const MuonServiceProxy * theProxyService
void TrackerSeedCleaner::setEvent ( const edm::Event event)
virtual

setEvent

Definition at line 53 of file TrackerSeedCleaner.cc.

Referenced by TSGFromL2Muon::produce().

54 {
55  event.getByToken(beamspotToken_, bsHandle_);
56 }
edm::Handle< reco::BeamSpot > bsHandle_
edm::EDGetTokenT< reco::BeamSpot > beamspotToken_

Member Data Documentation

edm::EDGetTokenT<reco::BeamSpot> TrackerSeedCleaner::beamspotToken_
private

Definition at line 66 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().

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

Definition at line 65 of file TrackerSeedCleaner.h.

std::string TrackerSeedCleaner::builderName_
private

Definition at line 69 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().

bool TrackerSeedCleaner::cleanBySharedHits
private

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().

edm::InputTag TrackerSeedCleaner::theBeamSpotTag
private

Definition at line 64 of file TrackerSeedCleaner.h.

Referenced by 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.

RedundantSeedCleaner* TrackerSeedCleaner::theRedundantCleaner
private

Definition at line 67 of file TrackerSeedCleaner.h.

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

Definition at line 70 of file TrackerSeedCleaner.h.

bool TrackerSeedCleaner::useDirection_Cleaner
private

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().

bool TrackerSeedCleaner::usePt_Cleaner
private

Definition at line 71 of file TrackerSeedCleaner.h.

Referenced by TrackerSeedCleaner().