Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:26:00

0001 /** \class MuonDetCleaner
0002  *
0003  * Clean collections of hits in muon detectors (CSC, DT and RPC)
0004  * for original Zmumu event and "embedded" simulated tau decay products
0005  * 
0006  * \author Christian Veelken, LLR
0007  *
0008  * 
0009  *
0010  * 
0011  *
0012  * Clean Up from STefan Wayand, KIT
0013  * 
0014  */
0015 #ifndef TauAnalysis_MCEmbeddingTools_MuonDetCleaner_H
0016 #define TauAnalysis_MCEmbeddingTools_MuonDetCleaner_H
0017 
0018 #include "FWCore/Framework/interface/stream/EDProducer.h"
0019 #include "FWCore/Framework/interface/Event.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "DataFormats/Common/interface/Handle.h"
0022 #include "FWCore/Framework/interface/MakerMacros.h"
0023 
0024 #include "DataFormats/Common/interface/RangeMap.h"
0025 #include "DataFormats/Common/interface/OwnVector.h"
0026 
0027 #include "DataFormats/PatCandidates/interface/Muon.h"
0028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0029 
0030 #include <string>
0031 #include <vector>
0032 #include <map>
0033 
0034 template <typename T1, typename T2>
0035 class MuonDetCleaner : public edm::stream::EDProducer<> {
0036 public:
0037   explicit MuonDetCleaner(const edm::ParameterSet&);
0038   ~MuonDetCleaner() override;
0039 
0040 private:
0041   void produce(edm::Event&, const edm::EventSetup&) override;
0042 
0043   typedef edm::RangeMap<T1, edm::OwnVector<T2> > RecHitCollection;
0044   void fillVetoHits(const TrackingRecHit&, std::vector<uint32_t>*);
0045 
0046   uint32_t getRawDetId(const T2&);
0047   bool checkrecHit(const TrackingRecHit&);
0048 
0049   const edm::EDGetTokenT<edm::View<pat::Muon> > mu_input_;
0050 
0051   std::map<std::string, edm::EDGetTokenT<RecHitCollection> > inputs_;
0052 };
0053 
0054 template <typename T1, typename T2>
0055 MuonDetCleaner<T1, T2>::MuonDetCleaner(const edm::ParameterSet& iConfig)
0056     : mu_input_(consumes<edm::View<pat::Muon> >(iConfig.getParameter<edm::InputTag>("MuonCollection"))) {
0057   std::vector<edm::InputTag> inCollections = iConfig.getParameter<std::vector<edm::InputTag> >("oldCollection");
0058   for (const auto& inCollection : inCollections) {
0059     inputs_[inCollection.instance()] = consumes<RecHitCollection>(inCollection);
0060     produces<RecHitCollection>(inCollection.instance());
0061   }
0062 }
0063 
0064 template <typename T1, typename T2>
0065 MuonDetCleaner<T1, T2>::~MuonDetCleaner() {
0066   // nothing to be done yet...
0067 }
0068 
0069 template <typename T1, typename T2>
0070 void MuonDetCleaner<T1, T2>::produce(edm::Event& iEvent, const edm::EventSetup& es) {
0071   std::map<T1, std::vector<T2> > recHits_output;  // This data format is easyer to handle
0072   std::vector<uint32_t> vetoHits;
0073 
0074   // First fill the veto RecHits colletion with the Hits from the input muons
0075   edm::Handle<edm::View<pat::Muon> > muonHandle;
0076   iEvent.getByToken(mu_input_, muonHandle);
0077   edm::View<pat::Muon> muons = *muonHandle;
0078   for (edm::View<pat::Muon>::const_iterator iMuon = muons.begin(); iMuon != muons.end(); ++iMuon) {
0079     const reco::Track* track = nullptr;
0080     if (iMuon->isGlobalMuon())
0081       track = iMuon->outerTrack().get();
0082     else if (iMuon->isStandAloneMuon())
0083       track = iMuon->outerTrack().get();
0084     else if (iMuon->isRPCMuon())
0085       track = iMuon->innerTrack().get();  // To add, try to access the rpc track
0086     else if (iMuon->isTrackerMuon())
0087       track = iMuon->innerTrack().get();
0088     else {
0089       edm::LogError("TauEmbedding") << "The imput muon: " << (*iMuon)
0090                                     << " must be either global or does or be tracker muon";
0091       assert(0);
0092     }
0093 
0094     for (trackingRecHit_iterator hitIt = track->recHitsBegin(); hitIt != track->recHitsEnd(); ++hitIt) {
0095       const TrackingRecHit& murechit = **hitIt;  // Base class for all rechits
0096       if (!(murechit).isValid())
0097         continue;
0098       if (!checkrecHit(murechit))
0099         continue;                         // Check if the hit belongs to a specifc detector section
0100       fillVetoHits(murechit, &vetoHits);  // Go back to the very basic rechits
0101     }
0102   }
0103 
0104   // Now this can also handle different instance
0105   for (auto input_ : inputs_) {
0106     // Second read in the RecHit Colltection which is to be replaced, without the vetoRecHits
0107     typedef edm::Handle<RecHitCollection> RecHitCollectionHandle;
0108     RecHitCollectionHandle RecHitinput;
0109     iEvent.getByToken(input_.second, RecHitinput);
0110     for (typename RecHitCollection::const_iterator recHit = RecHitinput->begin(); recHit != RecHitinput->end();
0111          ++recHit) {  // loop over the basic rec hit collection (DT CSC or RPC)
0112       //if (find(vetoHits.begin(),vetoHits.end(),getRawDetId(*recHit)) == vetoHits.end()) continue; // For the invertec selcetion
0113       if (find(vetoHits.begin(), vetoHits.end(), getRawDetId(*recHit)) != vetoHits.end())
0114         continue;  // If the hit is not in the
0115       T1 detId(getRawDetId(*recHit));
0116       recHits_output[detId].push_back(*recHit);
0117     }
0118 
0119     // Last step savet the output in the CMSSW Data Format
0120     std::unique_ptr<RecHitCollection> output(new RecHitCollection());
0121     for (typename std::map<T1, std::vector<T2> >::const_iterator recHit = recHits_output.begin();
0122          recHit != recHits_output.end();
0123          ++recHit) {
0124       output->put(recHit->first, recHit->second.begin(), recHit->second.end());
0125     }
0126     output->post_insert();
0127     iEvent.put(std::move(output), input_.first);
0128   }
0129 }
0130 
0131 template <typename T1, typename T2>
0132 void MuonDetCleaner<T1, T2>::fillVetoHits(const TrackingRecHit& rh, std::vector<uint32_t>* HitsList) {
0133   std::vector<const TrackingRecHit*> rh_components = rh.recHits();
0134   if (rh_components.empty()) {
0135     HitsList->push_back(rh.rawId());
0136   } else {
0137     for (std::vector<const TrackingRecHit*>::const_iterator rh_component = rh_components.begin();
0138          rh_component != rh_components.end();
0139          ++rh_component) {
0140       fillVetoHits(**rh_component, HitsList);
0141     }
0142   }
0143 }
0144 
0145 #endif