Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:25:31

0001 /** \class RPCDigiReader
0002  *  Analyse the RPC digitizer (derived from R. Bellan DTDigiReader. 
0003  *  
0004  *  \authors: M. Maggi -- INFN Bari
0005  */
0006 
0007 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 
0014 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0015 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0016 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
0017 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0018 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
0019 #include "DataFormats/RPCDigi/interface/RPCDigiCollection.h"
0020 #include "SimDataFormats/RPCDigiSimLink/interface/RPCDigiSimLink.h"
0021 #include "DataFormats/Common/interface/DetSetVector.h"
0022 #include <map>
0023 #include <set>
0024 
0025 #include "DataFormats/Common/interface/DetSet.h"
0026 
0027 #include <iostream>
0028 
0029 class RPCDigiReader : public edm::one::EDAnalyzer<> {
0030 public:
0031   explicit RPCDigiReader(const edm::ParameterSet& pset)
0032       : label_(pset.getUntrackedParameter<std::string>("label")),
0033         tokGeom_(esConsumes<RPCGeometry, MuonGeometryRecord>()),
0034         digiToken_(consumes<RPCDigiCollection>(edm::InputTag(label_))),
0035         hitsToken_(consumes<edm::PSimHitContainer>(edm::InputTag("g4SimHits", "MuonRPCHits"))),
0036         linkToken_(consumes<edm::DetSetVector<RPCDigiSimLink> >(edm::InputTag("muonRPCDigis", "RPCDigiSimLink"))) {}
0037 
0038   ~RPCDigiReader() override = default;
0039 
0040   void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) override {
0041     edm::LogVerbatim("RPCDump") << "--- Run: " << event.id().run() << " Event: " << event.id().event();
0042 
0043     const auto& rpcDigis = event.getHandle(digiToken_);
0044     const auto& simHits = event.getHandle(hitsToken_);
0045     const auto& thelinkDigis = event.getHandle(linkToken_);
0046 
0047     const auto& pDD = eventSetup.getHandle(tokGeom_);
0048 
0049     RPCDigiCollection::DigiRangeIterator detUnitIt;
0050     for (detUnitIt = rpcDigis->begin(); detUnitIt != rpcDigis->end(); ++detUnitIt) {
0051       const RPCDetId& id = (*detUnitIt).first;
0052       const RPCRoll* roll = dynamic_cast<const RPCRoll*>(pDD->roll(id));
0053       const RPCDigiCollection::Range& range = (*detUnitIt).second;
0054 
0055       //     if(id.rawId() != 637567293) continue;
0056 
0057       // RPCDetId print-out
0058       edm::LogVerbatim("RPCDump") << "--------------";
0059       edm::LogVerbatim("RPCDump") << "id: " << id.rawId() << " number of strip " << roll->nstrips();
0060 
0061       // Loop over the digis of this DetUnit
0062       for (RPCDigiCollection::const_iterator digiIt = range.first; digiIt != range.second; ++digiIt) {
0063         edm::LogVerbatim("RPCDump") << " digi " << *digiIt;
0064         if (digiIt->strip() < 1 || digiIt->strip() > roll->nstrips()) {
0065           edm::LogVerbatim("RPCDump") << " XXXXXXXXXXXXX Problemt with " << id;
0066         }
0067         for (std::vector<PSimHit>::const_iterator simHit = simHits->begin(); simHit != simHits->end(); simHit++) {
0068           RPCDetId rpcId((*simHit).detUnitId());
0069           if (rpcId == id && abs((*simHit).particleType()) == 13) {
0070             edm::LogVerbatim("RPCDump") << "entry: " << (*simHit).entryPoint() << "\nexit: " << (*simHit).exitPoint()
0071                                         << "\nTOF: " << (*simHit).timeOfFlight();
0072           }
0073         }
0074       }  // for digis in layer
0075     }    // for layers
0076 
0077     for (edm::DetSetVector<RPCDigiSimLink>::const_iterator itlink = thelinkDigis->begin();
0078          itlink != thelinkDigis->end();
0079          itlink++) {
0080       for (edm::DetSet<RPCDigiSimLink>::const_iterator digi_iter = itlink->data.begin();
0081            digi_iter != itlink->data.end();
0082            ++digi_iter) {
0083         int ev = digi_iter->getEventId().event();
0084         int detid = digi_iter->getDetUnitId();
0085         float xpos = digi_iter->getEntryPoint().x();
0086         int strip = digi_iter->getStrip();
0087         int bx = digi_iter->getBx();
0088 
0089         edm::LogVerbatim("RPCDump") << "DetUnit: " << detid << "  "
0090                                     << "Event ID: " << ev << "  "
0091                                     << "Pos X: " << xpos << "  "
0092                                     << "Strip: " << strip << "  "
0093                                     << "Bx: " << bx;
0094       }
0095     }
0096 
0097     edm::LogVerbatim("RPCDump") << "--------------";
0098   }
0099 
0100 private:
0101   const std::string label_;
0102   const edm::ESGetToken<RPCGeometry, MuonGeometryRecord> tokGeom_;
0103   const edm::EDGetTokenT<RPCDigiCollection> digiToken_;
0104   const edm::EDGetTokenT<edm::PSimHitContainer> hitsToken_;
0105   const edm::EDGetTokenT<edm::DetSetVector<RPCDigiSimLink> > linkToken_;
0106 };
0107 
0108 #include "FWCore/Framework/interface/MakerMacros.h"
0109 DEFINE_FWK_MODULE(RPCDigiReader);