CMS 3D CMS Logo

TrackerCleaner.h
Go to the documentation of this file.
1 
13 #ifndef TauAnalysis_MCEmbeddingTools_TrackerCleaner_H
14 #define TauAnalysis_MCEmbeddingTools_TrackerCleaner_H
15 
23 
27 
32 
33 #include <string>
34 #include <iostream>
35 #include <map>
36 
37 template <typename T>
39 public:
40  explicit TrackerCleaner(const edm::ParameterSet&);
41  ~TrackerCleaner() override;
42 
43 private:
44  void produce(edm::Event&, const edm::EventSetup&) override;
45 
48 
49  std::map<std::string, edm::EDGetTokenT<TrackClusterCollection> > inputs_;
50 
51  bool match_rechit_type(const TrackingRecHit& murechit);
52 };
53 
54 template <typename T>
56  : mu_input_(consumes<edm::View<pat::Muon> >(iConfig.getParameter<edm::InputTag>("MuonCollection")))
57 
58 {
59  std::vector<edm::InputTag> inCollections = iConfig.getParameter<std::vector<edm::InputTag> >("oldCollection");
60  for (const auto& inCollection : inCollections) {
61  inputs_[inCollection.instance()] = consumes<TrackClusterCollection>(inCollection);
62  produces<TrackClusterCollection>(inCollection.instance());
63  }
64 }
65 
66 template <typename T>
68  // nothing to be done yet...
69 }
70 
71 template <typename T>
73  using namespace edm;
74 
76  iEvent.getByToken(mu_input_, muonHandle);
77  edm::View<pat::Muon> muons = *muonHandle;
78 
79  for (auto input_ : inputs_) {
81  iEvent.getByToken(input_.second, inputClusters);
82 
83  std::vector<bool> vetodClusters;
84 
85  vetodClusters.resize(inputClusters->dataSize(), false);
86 
87  for (edm::View<pat::Muon>::const_iterator iMuon = muons.begin(); iMuon != muons.end(); ++iMuon) {
88  if (!iMuon->isGlobalMuon())
89  continue;
90  const reco::Track* mutrack = iMuon->globalTrack().get();
91  // reco::Track *mutrack = new reco::Track(*(iMuon->innerTrack() ));
92  for (trackingRecHit_iterator hitIt = mutrack->recHitsBegin(); hitIt != mutrack->recHitsEnd(); ++hitIt) {
93  const TrackingRecHit& murechit = **hitIt;
94  if (!(murechit).isValid())
95  continue;
96 
97  if (match_rechit_type(murechit)) {
98  auto& thit = reinterpret_cast<BaseTrackerRecHit const&>(murechit);
99  auto const& cluster = thit.firstClusterRef();
100  vetodClusters[cluster.key()] = true;
101  }
102  }
103  }
104  std::unique_ptr<TrackClusterCollection> output(new TrackClusterCollection());
105 
106  int idx = 0;
107  for (typename TrackClusterCollection::const_iterator clustSet = inputClusters->begin();
108  clustSet != inputClusters->end();
109  ++clustSet) {
110  DetId detIdObject(clustSet->detId());
111  typename TrackClusterCollection::FastFiller spc(*output, detIdObject);
112  for (typename edmNew::DetSet<T>::const_iterator clustIt = clustSet->begin(); clustIt != clustSet->end();
113  ++clustIt) {
114  idx++;
115  if (vetodClusters[idx - 1])
116  continue;
117  //if (!vetodClusters[idx-1]) continue; for inverted selction
118  spc.push_back(*clustIt);
119  }
120  }
121  iEvent.put(std::move(output), input_.first);
122  }
123 }
124 #endif
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::EDGetTokenT< edm::View< pat::Muon > > mu_input_
data_type const * const_iterator
Definition: DetSetNew.h:31
TrackerCleaner(const edm::ParameterSet &)
void produce(edm::Event &, const edm::EventSetup &) override
~TrackerCleaner() override
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
Definition: Muon.py:1
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:91
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:88
std::map< std::string, edm::EDGetTokenT< TrackClusterCollection > > inputs_
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
bool match_rechit_type(const TrackingRecHit &murechit)
def move(src, dest)
Definition: eostools.py:511
edmNew::DetSetVector< T > TrackClusterCollection