CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastTSGFromIOHit.cc
Go to the documentation of this file.
2 
17 
18 
20 {
21  theCategory = "FastSimulation|Muons||FastTSGFromIOHit";
22 
23  thePtCut = iConfig.getParameter<double>("PtCut");
24 
25  theSeedCollectionLabels = iConfig.getParameter<std::vector<edm::InputTag> >("SeedCollectionLabels");
26 
27  theSimTrackCollectionLabel = iConfig.getParameter<edm::InputTag>("SimTrackCollectionLabel");
28 
29 }
30 
32 {
33 
34  LogTrace(theCategory) << " FastTSGFromIOHit dtor called ";
35 
36 }
37 
38 void FastTSGFromIOHit::trackerSeeds(const TrackCand& staMuon, const TrackingRegion& region, std::vector<TrajectorySeed> & result) {
39 
40  // Retrieve the Monte Carlo truth (SimTracks)
43 
44  // Retrieve Seed collection
45  unsigned seedCollections = theSeedCollectionLabels.size();
46  std::vector<edm::Handle<edm::View<TrajectorySeed> > > theSeeds;
47  theSeeds.resize(seedCollections);
48  unsigned seed_size = 0;
49  for ( unsigned iseed=0; iseed<seedCollections; ++iseed ) {
50  getEvent()->getByLabel(theSeedCollectionLabels[iseed], theSeeds[iseed]);
51  seed_size += theSeeds[iseed]->size();
52  }
53 
54  // Make a ref to l2 muon
55  reco::TrackRef muRef(staMuon.second);
56 
57  // Cut on muons with low momenta
58  if ( muRef->pt() < thePtCut
59  || muRef->innerMomentum().Rho() < thePtCut
60  || muRef->innerMomentum().R() < 2.5 ){
61  }return;
62 
63  // Copy the collection of seeds (ahem, this is time consuming!)
64  std::vector<TrajectorySeed> tkSeeds;
65  std::set<unsigned> tkIds;
66  tkSeeds.reserve(seed_size);
67  for ( unsigned iseed=0; iseed<seedCollections; ++iseed ) {
68  edm::Handle<edm::View<TrajectorySeed> > aSeedCollection = theSeeds[iseed];
69  unsigned nSeeds = aSeedCollection->size();
70  for (unsigned seednr = 0; seednr < nSeeds; ++seednr) {
71 
72  // The seed
73  const BasicTrajectorySeed* aSeed = &((*aSeedCollection)[seednr]);
74 
75  // Find the first hit of the Seed
76  TrajectorySeed::range theSeedingRecHitRange = aSeed->recHits();
77  const SiTrackerGSMatchedRecHit2D * theFirstSeedingRecHit =
78  (const SiTrackerGSMatchedRecHit2D*) (&(*(theSeedingRecHitRange.first)));
79 
80  // The SimTrack id associated to that recHit
81  int simTrackId = theFirstSeedingRecHit->simtrackId();
82 
83  // Track already associated to a seed
84  std::set<unsigned>::iterator tkId = tkIds.find(simTrackId);
85  if( tkId != tkIds.end() ) continue;
86 
87  const SimTrack& theSimTrack = (*theSimTracks)[simTrackId];
88 
89  const RectangularEtaPhiTrackingRegion & regionRef = dynamic_cast<const RectangularEtaPhiTrackingRegion & > (region);
90 
91  if( clean(muRef,regionRef,aSeed,theSimTrack) ) tkSeeds.push_back(*aSeed);
92  tkIds.insert(simTrackId);
93 
94  } // End loop on seeds
95 
96  } // End loop on seed collections
97 
98  // Now create the Muon Trajectory Seed
99  unsigned int is=0;
100  unsigned int isMax=tkSeeds.size();
101  for (;is!=isMax;++is){
102  result.push_back( L3MuonTrajectorySeed(tkSeeds[is], muRef));
103  } // End of tk seed loop
104 
105 }
106 
107 bool
110  const BasicTrajectorySeed* aSeed,
111  const SimTrack& theSimTrack)
112 {
113  // return true;
114  //}
115 
116  // Eta cleaner
117  const PixelRecoRange<float>& etaRange = region.etaRange() ;
118  // return true;
119  //}
120 
121  double etaSeed = theSimTrack.momentum().Eta();
122  double etaLimit = (fabs(fabs(etaRange.max())-fabs(etaRange.mean())) <0.05) ?
123  0.05 : fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) ;
124  bool inEtaRange =
125  etaSeed >= (etaRange.mean() - etaLimit) &&
126  etaSeed <= (etaRange.mean() + etaLimit) ;
127  if ( !inEtaRange ) return false;
128 
129  // Phi cleaner
130  const TkTrackingRegionsMargin<float>& phiMargin = region.phiMargin();
131  double phiSeed = theSimTrack.momentum().Phi();
132  double phiLimit = (phiMargin.right() < 0.05 ) ? 0.05 : phiMargin.right();
133  bool inPhiRange =
134  (fabs(deltaPhi(phiSeed,double(region.direction().phi()))) < phiLimit );
135  if ( !inPhiRange ) return false;
136 
137  // pt cleaner
138  double ptSeed = std::sqrt(theSimTrack.momentum().Perp2());
139  double ptMin = (region.ptMin()>3.5) ? 3.5: region.ptMin();
140  bool inPtRange = ptSeed >= ptMin && ptSeed<= 2*(muRef->pt());
141  return inPtRange;
142 
143 }
T getParameter(std::string const &) const
T max() const
FastTSGFromIOHit(const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
constructor
std::pair< const Trajectory *, reco::TrackRef > TrackCand
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
edm::InputTag theSimTrackCollectionLabel
virtual ~FastTSGFromIOHit()
destructor
std::vector< edm::InputTag > theSeedCollectionLabels
GlobalVector const & direction() const
the direction around which region is constructed
void trackerSeeds(const TrackCand &, const TrackingRegion &, std::vector< TrajectorySeed > &)
generate seed(s) for a track
T sqrt(T t)
Definition: SSEVec.h:48
const edm::Event * getEvent() const
tuple result
Definition: query.py:137
std::pair< const_iterator, const_iterator > range
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:413
#define LogTrace(id)
int iseed
Definition: AMPTWrapper.h:124
T mean() const
range recHits() const
float ptMin() const
minimal pt of interest
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
bool clean(reco::TrackRef muRef, const RectangularEtaPhiTrackingRegion &region, const BasicTrajectorySeed *aSeed, const SimTrack &theSimTrack)
std::string theCategory
const Range & etaRange() const
allowed eta range [eta_min, eta_max] interval