CMS 3D CMS Logo

FastTSGFromL2Muon.cc
Go to the documentation of this file.
2 
5 
10 
13 
14 #include <set>
15 
17  : thePtCut(cfg.getParameter<double>("PtCut")),
18  theL2CollectionLabel(cfg.getParameter<edm::InputTag>("MuonCollectionLabel")),
19  theSeedCollectionLabels(cfg.getParameter<std::vector<edm::InputTag> >("SeedCollectionLabels")),
20  theSimTrackCollectionLabel(cfg.getParameter<edm::InputTag>("SimTrackCollectionLabel")),
21  simTrackToken_(consumes<edm::SimTrackContainer>(theSimTrackCollectionLabel)),
22  l2TrackToken_(consumes<reco::TrackCollection>(theL2CollectionLabel)) {
23  produces<L3MuonTrajectorySeedCollection>();
24 
25  for (auto& seed : theSeedCollectionLabels)
26  seedToken_.emplace_back(consumes<edm::View<TrajectorySeed> >(seed));
27 
28  edm::ParameterSet regionBuilderPSet = cfg.getParameter<edm::ParameterSet>("MuonTrackingRegionBuilder");
29  theRegionBuilder = std::make_unique<MuonTrackingRegionBuilder>(regionBuilderPSet, consumesCollector());
30 }
31 
33  //region builder
34 }
35 
37  // Initialize the output product
38  std::unique_ptr<L3MuonTrajectorySeedCollection> result(new L3MuonTrajectorySeedCollection());
39 
40  // Region builder
41  theRegionBuilder->setEvent(ev, es);
42 
43  // Retrieve the Monte Carlo truth (SimTracks)
44  const edm::Handle<edm::SimTrackContainer>& theSimTracks = ev.getHandle(simTrackToken_);
45 
46  // Retrieve L2 muon collection
47  const edm::Handle<reco::TrackCollection>& l2muonH = ev.getHandle(l2TrackToken_);
48 
49  // Retrieve Seed collection
50  unsigned seedCollections = theSeedCollectionLabels.size();
51  std::vector<edm::Handle<edm::View<TrajectorySeed> > > theSeeds;
52  theSeeds.resize(seedCollections);
53  unsigned seed_size = 0;
54  for (unsigned iseed = 0; iseed < seedCollections; ++iseed) {
55  ev.getByToken(seedToken_[iseed], theSeeds[iseed]);
56  seed_size += theSeeds[iseed]->size();
57  }
58 
59  // Loop on L2 muons
60  unsigned int imu = 0;
61  unsigned int imuMax = l2muonH->size();
62  edm::LogVerbatim("FastTSGFromL2Muon") << "Found " << imuMax << " L2 muons";
63  for (; imu != imuMax; ++imu) {
64  // Make a ref to l2 muon
65  reco::TrackRef muRef(l2muonH, imu);
66 
67  // Cut on muons with low momenta
68  if (muRef->pt() < thePtCut || muRef->innerMomentum().Rho() < thePtCut || muRef->innerMomentum().R() < 2.5)
69  continue;
70 
71  // Define the region of interest
72  std::unique_ptr<RectangularEtaPhiTrackingRegion> region = theRegionBuilder->region(muRef);
73 
74  // Copy the collection of seeds (ahem, this is time consuming!)
75  std::vector<TrajectorySeed> tkSeeds;
76  std::set<unsigned> tkIds;
77  tkSeeds.reserve(seed_size);
78  for (unsigned iseed = 0; iseed < seedCollections; ++iseed) {
79  edm::Handle<edm::View<TrajectorySeed> > aSeedCollection = theSeeds[iseed];
80  unsigned nSeeds = aSeedCollection->size();
81  for (unsigned seednr = 0; seednr < nSeeds; ++seednr) {
82  // The seed
83  const BasicTrajectorySeed* aSeed = &((*aSeedCollection)[seednr]);
84 
85  // The SimTrack id associated to the first hit of the Seed
86  int simTrackId = static_cast<FastTrackerRecHit const&>(*aSeed->recHits().begin()).simTrackId(0);
87 
88  // Track already associated to a seed
89  std::set<unsigned>::iterator tkId = tkIds.find(simTrackId);
90  if (tkId != tkIds.end())
91  continue;
92 
93  const SimTrack& theSimTrack = (*theSimTracks)[simTrackId];
94 
95  if (clean(muRef, region.get(), aSeed, theSimTrack))
96  tkSeeds.push_back(*aSeed);
97  tkIds.insert(simTrackId);
98 
99  } // End loop on seeds
100 
101  } // End loop on seed collections
102 
103  // Now create the Muon Trajectory Seed
104  unsigned int is = 0;
105  unsigned int isMax = tkSeeds.size();
106  for (; is != isMax; ++is) {
107  result->push_back(L3MuonTrajectorySeed(tkSeeds[is], muRef));
108  } // End of tk seed loop
109 
110  } // End of l2 muon loop
111 
112  edm::LogVerbatim("FastTSGFromL2Muon") << "Found " << result->size() << " seeds for muons";
113 
114  //put in the event
115  ev.put(std::move(result));
116 }
117 
120  const BasicTrajectorySeed* aSeed,
121  const SimTrack& theSimTrack) {
122  // Eta cleaner
123  const PixelRecoRange<float>& etaRange = region->etaRange();
124  double etaSeed = theSimTrack.momentum().Eta();
125  double etaLimit = (fabs(fabs(etaRange.max()) - fabs(etaRange.mean())) < 0.05)
126  ? 0.05
127  : fabs(fabs(etaRange.max()) - fabs(etaRange.mean()));
128  bool inEtaRange = etaSeed >= (etaRange.mean() - etaLimit) && etaSeed <= (etaRange.mean() + etaLimit);
129  if (!inEtaRange)
130  return false;
131 
132  // Phi cleaner
133  const TkTrackingRegionsMargin<float>& phiMargin = region->phiMargin();
134  double phiSeed = theSimTrack.momentum().Phi();
135  double phiLimit = (phiMargin.right() < 0.05) ? 0.05 : phiMargin.right();
136  bool inPhiRange = (fabs(deltaPhi(phiSeed, double(region->direction().phi()))) < phiLimit);
137  if (!inPhiRange)
138  return false;
139 
140  // pt cleaner
141  double ptSeed = std::sqrt(theSimTrack.momentum().Perp2());
142  double ptMin = (region->ptMin() > 3.5) ? 3.5 : region->ptMin();
143  bool inPtRange = ptSeed >= ptMin && ptSeed <= 2 * (muRef->pt());
144  return inPtRange;
145 }
Log< level::Info, true > LogVerbatim
const double thePtCut
const edm::EDGetTokenT< edm::SimTrackContainer > simTrackToken_
void produce(edm::Event &ev, const edm::EventSetup &es) override
T begin() const
Definition: Range.h:15
RecHitRange recHits() const
std::vector< L3MuonTrajectorySeed > L3MuonTrajectorySeedCollection
constexpr float ptMin
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< edm::EDGetTokenT< edm::View< TrajectorySeed > > > seedToken_
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
const edm::EDGetTokenT< reco::TrackCollection > l2TrackToken_
T sqrt(T t)
Definition: SSEVec.h:19
std::unique_ptr< MuonTrackingRegionBuilder > theRegionBuilder
bool clean(reco::TrackRef muRef, RectangularEtaPhiTrackingRegion *region, const BasicTrajectorySeed *aSeed, const SimTrack &theSimTrack)
FastTSGFromL2Muon(const edm::ParameterSet &cfg)
int iseed
Definition: AMPTWrapper.h:134
const std::vector< edm::InputTag > theSeedCollectionLabels
fixed size matrix
HLT enums.
void beginRun(edm::Run const &run, edm::EventSetup const &es) override
std::vector< SimTrack > SimTrackContainer
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45