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)
 constructor More...
 
virtual ~TrackerSeedCleaner ()
 destructor More...
 

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(), AlCaHLTBitMon_QueryRunRegistry::string, theBeamSpotTag, useDirection_Cleaner, and usePt_Cleaner.

41  : theProxyService(0),theEvent(0) {
42  builderName_ = pset.getParameter<std::string>("TTRHBuilder");
43  theBeamSpotTag = pset.getParameter<edm::InputTag>("beamSpot");
44  useDirection_Cleaner = pset.getParameter<bool>("directionCleaner");
45  usePt_Cleaner = pset.getParameter<bool>("ptCleaner");
46  cleanBySharedHits = pset.getParameter<bool>("cleanerFromSharedHits");
47  }
T getParameter(std::string const &) const
edm::InputTag theBeamSpotTag
const edm::Event * theEvent
const MuonServiceProxy * theProxyService
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 63 of file TrackerSeedCleaner.cc.

References SiPixelRawToDigiRegional_cfi::deltaPhi, PV3DBase< T, PVType, FrameType >::eta(), MuonErrorMatrixValues_cff::etaRange, TrajectoryStateOnSurface::freeState(), TrajectoryStateOnSurface::globalMomentum(), LogDebug, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), PtMinSelector_cfg::ptMin, query::result, evf::utils::state, 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().

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

intizialization

Definition at line 45 of file TrackerSeedCleaner.cc.

Referenced by TSGFromL2Muon::beginRun().

45  {
46 
47  theProxyService = service;
48 
50 }
RedundantSeedCleaner * theRedundantCleaner
const MuonServiceProxy * theProxyService
void TrackerSeedCleaner::setEvent ( const edm::Event event)
virtual

setEvent

Definition at line 55 of file TrackerSeedCleaner.cc.

Referenced by TSGFromL2Muon::produce().

56 {
57  event.getByLabel(theBeamSpotTag, bsHandle_);
58 }
edm::InputTag theBeamSpotTag
edm::Handle< reco::BeamSpot > bsHandle_

Member Data Documentation

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().