CMS 3D CMS Logo

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