CMS 3D CMS Logo

TrackerSeedCleaner.cc
Go to the documentation of this file.
1 /*
2  * \class TrackerSeedCleaner
3  * Reference class for seeds cleaning
4  * Seeds Cleaner based on sharedHits cleaning, direction cleaning and pt cleaning
5  \author A. Grelli - Purdue University, Pavia University
6  */
7 
9 
10 //---------------
11 // C++ Headers --
12 //---------------
13 #include <vector>
14 
15 //-------------------------------
16 // Collaborating Class Headers --
17 //-------------------------------
20 
24 
32 
36 
37 using namespace std;
38 using namespace edm;
39 
40 //
41 // inizialization
42 //
44 
45  theProxyService = service;
46 
47  theRedundantCleaner = new RedundantSeedCleaner();
48 }
49 
50 //
51 //
52 //
54 {
55  event.getByToken(beamspotToken_, bsHandle_);
56 }
57 
58 //
59 // clean seeds
60 //
62 
63 
64  // call the shared input cleaner
65  if(cleanBySharedHits) theRedundantCleaner->define(seeds);
66 
67  theProxyService->eventSetup().get<TransientRecHitRecord>().get(builderName_,theTTRHBuilder);
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.empty() && (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 }
159 
#define LogDebug(id)
double z0() const
z coordinate
Definition: BeamSpot.h:68
T perp() const
Definition: PV3DBase.h:72
PixelRecoRange< float > Range
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
virtual void setEvent(const edm::Event &)
setEvent
GlobalVector const & direction() const
the direction around which region is constructed
FreeTrajectoryState const * freeState(bool withErrors=true) const
std::vector< TrajectorySeed > tkSeeds
T z() const
Definition: PV3DBase.h:64
std::shared_ptr< TrackingRecHit const > RecHitPointer
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
float ptMin() const
minimal pt of interest
T eta() const
Definition: PV3DBase.h:76
HLT enums.
GlobalVector globalMomentum() const
double y0() const
y coordinate
Definition: BeamSpot.h:66
virtual void init(const MuonServiceProxy *service)
intizialization
T x() const
Definition: PV3DBase.h:62
Definition: event.py:1
virtual void clean(const reco::TrackRef &, const RectangularEtaPhiTrackingRegion &region, tkSeeds &)
the cleaner
double x0() const
x coordinate
Definition: BeamSpot.h:64
const Range & etaRange() const
allowed eta range [eta_min, eta_max] interval