CMS 3D CMS Logo

MuonDetCleaner.h
Go to the documentation of this file.
1 
15 #ifndef TauAnalysis_MCEmbeddingTools_MuonDetCleaner_H
16 #define TauAnalysis_MCEmbeddingTools_MuonDetCleaner_H
17 
18 
24 
27 
30 
31 #include <string>
32 #include <vector>
33 #include <map>
34 
35 template <typename T1, typename T2>
37 {
38  public:
39  explicit MuonDetCleaner(const edm::ParameterSet&);
40  ~MuonDetCleaner() override;
41 
42  private:
43  void produce(edm::Event&, const edm::EventSetup&) override;
44 
46  void fillVetoHits(const TrackingRecHit& , std::vector<uint32_t>* );
47 
48  uint32_t getRawDetId(const T2&);
49  bool checkrecHit(const TrackingRecHit&);
50 
52 
53  std::map<std::string, edm::EDGetTokenT<RecHitCollection > > inputs_;
54 
55 };
56 
57 template <typename T1, typename T2>
59  mu_input_(consumes<edm::View<pat::Muon> >(iConfig.getParameter<edm::InputTag>("MuonCollection")))
60 {
61  std::vector<edm::InputTag> inCollections = iConfig.getParameter<std::vector<edm::InputTag> >("oldCollection");
62  for (auto inCollection : inCollections){
63  inputs_[inCollection.instance()] = consumes<RecHitCollection >(inCollection);
64  produces<RecHitCollection>(inCollection.instance());
65  }
66 
67 }
68 
69 template <typename T1, typename T2>
71 {
72 // nothing to be done yet...
73 }
74 
75 template <typename T1, typename T2>
77 {
78  std::map<T1, std::vector<T2> > recHits_output; // This data format is easyer to handle
79  std::vector<uint32_t> vetoHits;
80 
81  // First fill the veto RecHits colletion with the Hits from the input muons
83  iEvent.getByToken(mu_input_, muonHandle);
84  edm::View<pat::Muon> muons = *muonHandle;
85  for (edm::View<pat::Muon>::const_iterator iMuon = muons.begin(); iMuon != muons.end(); ++iMuon) {
86  const reco::Track* track = nullptr;
87  if( iMuon->isGlobalMuon() ) track = iMuon->outerTrack().get();
88  else if ( iMuon->isStandAloneMuon() ) track = iMuon->outerTrack().get();
89  else if ( iMuon->isRPCMuon() ) track = iMuon->innerTrack().get(); // To add, try to access the rpc track
90  else if ( iMuon->isTrackerMuon() ) track = iMuon->innerTrack().get();
91  else {
92  edm::LogError("TauEmbedding")<<"The imput muon: "<<(*iMuon)<<" must be either global or does or be tracker muon";
93  assert(0);
94  }
95 
96 
97  for (trackingRecHit_iterator hitIt = track->recHitsBegin(); hitIt != track->recHitsEnd(); ++hitIt) {
98  const TrackingRecHit &murechit = **hitIt; // Base class for all rechits
99  if(!(murechit).isValid()) continue;
100  if (!checkrecHit(murechit)) continue; // Check if the hit belongs to a specifc detector section
101  fillVetoHits(murechit,&vetoHits); // Go back to the very basic rechits
102  }
103  }
104 
105 
106  // Now this can also handle different instance
107  for (auto input_ : inputs_){
108  // Second read in the RecHit Colltection which is to be replaced, without the vetoRecHits
109  typedef edm::Handle<RecHitCollection> RecHitCollectionHandle;
110  RecHitCollectionHandle RecHitinput;
111  iEvent.getByToken(input_.second, RecHitinput);
112  for ( typename RecHitCollection::const_iterator recHit = RecHitinput->begin(); recHit != RecHitinput->end(); ++recHit ) { // loop over the basic rec hit collection (DT CSC or RPC)
113  //if (find(vetoHits.begin(),vetoHits.end(),getRawDetId(*recHit)) == vetoHits.end()) continue; // For the invertec selcetion
114  if (find(vetoHits.begin(),vetoHits.end(),getRawDetId(*recHit)) != vetoHits.end()) continue; // If the hit is not in the
115  T1 detId(getRawDetId(*recHit));
116  recHits_output[detId].push_back(*recHit);
117  }
118 
119 
120 
121  // Last step savet the output in the CMSSW Data Format
122  std::unique_ptr<RecHitCollection> output(new RecHitCollection());
123  for ( typename std::map<T1, std::vector<T2> >::const_iterator recHit = recHits_output.begin(); recHit != recHits_output.end(); ++recHit ) {
124  output->put(recHit->first, recHit->second.begin(), recHit->second.end());
125  }
126  output->post_insert();
127  iEvent.put(std::move(output),input_.first);
128  }
129 }
130 
131 
132 template <typename T1, typename T2>
133 void MuonDetCleaner<T1,T2>::fillVetoHits(const TrackingRecHit& rh, std::vector<uint32_t>* HitsList)
134 {
135  std::vector<const TrackingRecHit*> rh_components = rh.recHits();
136  if ( rh_components.empty() ) {
137  HitsList->push_back(rh.rawId());
138  }
139  else {
140  for ( std::vector<const TrackingRecHit*>::const_iterator rh_component = rh_components.begin(); rh_component != rh_components.end(); ++rh_component ) {
141  fillVetoHits(**rh_component, HitsList);
142  }
143  }
144 }
145 
146 
147 
148 
149 #endif
150 
151 
152 
153 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
uint32_t getRawDetId(const T2 &)
bool checkrecHit(const TrackingRecHit &)
edm::RangeMap< T1, edm::OwnVector< T2 > > RecHitCollection
~MuonDetCleaner() override
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
Definition: HeavyIon.h:7
virtual std::vector< const TrackingRecHit * > recHits() const =0
Access to component RecHits (if any)
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int iEvent
Definition: GenABIO.cc:230
const_iterator begin() const
Definition: Muon.py:1
std::map< std::string, edm::EDGetTokenT< RecHitCollection > > inputs_
void produce(edm::Event &, const edm::EventSetup &) override
const edm::EDGetTokenT< edm::View< pat::Muon > > mu_input_
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
MuonDetCleaner(const edm::ParameterSet &)
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
void fillVetoHits(const TrackingRecHit &, std::vector< uint32_t > *)
const_iterator end() const
def move(src, dest)
Definition: eostools.py:510
id_type rawId() const
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109