Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:34

0001 #include <string>
0002 #include <iostream>
0003 #include <map>
0004 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/Framework/interface/Frameworkfwd.h"
0008 #include "DataFormats/SiPixelDetId/interface/PixelFEDChannel.h"
0009 #include "CalibTracker/Records/interface/SiPixelFEDChannelContainerESProducerRcd.h"
0010 
0011 class PixelFEDChannelCollectionMapTestReader : public edm::one::EDAnalyzer<> {
0012 public:
0013   typedef std::unordered_map<std::string, PixelFEDChannelCollection> PixelFEDChannelCollectionMap;
0014   explicit PixelFEDChannelCollectionMapTestReader(edm::ParameterSet const& p);
0015   ~PixelFEDChannelCollectionMapTestReader();
0016 
0017   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0018 
0019 private:
0020   virtual void analyze(const edm::Event& e, const edm::EventSetup& c) override;
0021 
0022   // ----------member data ---------------------------
0023   const bool printdebug_;
0024   const std::string formatedOutput_;
0025   edm::ESGetToken<PixelFEDChannelCollectionMap, SiPixelFEDChannelContainerESProducerRcd>
0026       pixelFEDChannelCollectionMapToken_;
0027 };
0028 
0029 PixelFEDChannelCollectionMapTestReader::PixelFEDChannelCollectionMapTestReader(edm::ParameterSet const& p)
0030     : printdebug_(p.getUntrackedParameter<bool>("printDebug", true)),
0031       formatedOutput_(p.getUntrackedParameter<std::string>("outputFile", "")),
0032       pixelFEDChannelCollectionMapToken_(
0033           esConsumes<PixelFEDChannelCollectionMap, SiPixelFEDChannelContainerESProducerRcd>()) {
0034   edm::LogInfo("PixelFEDChannelCollectionMapTestReader") << "PixelFEDChannelCollectionMapTestReader" << std::endl;
0035 }
0036 
0037 PixelFEDChannelCollectionMapTestReader::~PixelFEDChannelCollectionMapTestReader() {
0038   edm::LogInfo("PixelFEDChannelCollectionMapTestReader") << "~PixelFEDChannelCollectionMapTestReader " << std::endl;
0039 }
0040 
0041 void PixelFEDChannelCollectionMapTestReader::analyze(const edm::Event& e, const edm::EventSetup& context) {
0042   edm::LogInfo("PixelFEDChannelCollectionMapTestReader")
0043       << "### PixelFEDChannelCollectionMapTestReader::analyze  ###" << std::endl;
0044   edm::LogInfo("PixelFEDChannelCollectionMapTestReader") << " I AM IN RUN NUMBER " << e.id().run() << std::endl;
0045   edm::LogInfo("PixelFEDChannelCollectionMapTestReader") << " ---EVENT NUMBER " << e.id().event() << std::endl;
0046 
0047   edm::eventsetup::EventSetupRecordKey recordKey(
0048       edm::eventsetup::EventSetupRecordKey::TypeTag::findType("'SiPixelFEDChannelContainerESProducerRcd"));
0049 
0050   if (recordKey.type() == edm::eventsetup::EventSetupRecordKey::TypeTag()) {
0051     //record not found
0052     edm::LogInfo("PixelFEDChannelCollectionMapTestReader") << "Record \"SiPixelFEDChannelContainerESProducerRcd"
0053                                                            << "\" does not exist " << std::endl;
0054   }
0055 
0056   //this part gets the handle of the event source and the record (i.e. the Database)
0057   edm::ESHandle<PixelFEDChannelCollectionMap> PixelFEDChannelCollectionMapHandle =
0058       context.getHandle(pixelFEDChannelCollectionMapToken_);
0059   edm::LogInfo("PixelFEDChannelCollectionMapTestReader") << "got eshandle and context" << std::endl;
0060 
0061   const PixelFEDChannelCollectionMap* thePixelFEDChannelCollectionMap = PixelFEDChannelCollectionMapHandle.product();
0062   edm::LogInfo("PixelFEDChannelCollectionMapTestReader") << "got SiPixelFEDChannelContainer* " << std::endl;
0063   edm::LogInfo("PixelFEDChannelCollectionMapTestReader")
0064       << "print  pointer address : " << thePixelFEDChannelCollectionMap << std::endl;
0065   edm::LogInfo("PixelFEDChannelCollectionMapTestReader")
0066       << "Size: " << thePixelFEDChannelCollectionMap->size() << std::endl;
0067   edm::LogInfo("PixelFEDChannelCollectionMapTestReader") << "Content of my PixelFEDChanneCollectionlMap: " << std::endl;
0068 
0069   FILE* pFile = NULL;
0070   if (formatedOutput_ != "")
0071     pFile = fopen(formatedOutput_.c_str(), "w");
0072   if (pFile) {
0073     fprintf(pFile, "PixelFEDChannelCollectionMap::printAll() \n");
0074     fprintf(pFile,
0075             " ========================================================================================================="
0076             "========== \n");
0077 
0078     for (auto it = thePixelFEDChannelCollectionMap->begin(); it != thePixelFEDChannelCollectionMap->end(); ++it) {
0079       fprintf(pFile,
0080               " ======================================================================================================="
0081               "============ \n");
0082       fprintf(pFile, "run : %s \n ", (it->first).c_str());
0083       const auto& thePixelFEDChannels = it->second;
0084 
0085       //std:: cout << thePixelFEDChannels.size() << std::endl;
0086       // for (edmNew::DetSetVector<PixelFEDChannel>::const_iterator DSViter = thePixelFEDChannels.begin(); DSViter != thePixelFEDChannels.end(); DSViter++) {
0087       //       unsigned int theDetId = DSViter->id();
0088       //       std:: cout << theDetId << std::endl;
0089       // }
0090 
0091       for (const auto& disabledChannels : thePixelFEDChannels) {
0092         //loop over different PixelFED in a PixelFED vector (module)
0093         fprintf(pFile, "DetId : %i \n", disabledChannels.detId());
0094         for (const auto& ch : disabledChannels) {
0095           fprintf(pFile,
0096                   "fed : %i | link : %2i | roc_first : %2i | roc_last: %2i \n",
0097                   ch.fed,
0098                   ch.link,
0099                   ch.roc_first,
0100                   ch.roc_last);
0101           //std::cout <<  disabledChannels.detId() << " "<< ch.fed << " " << ch.link << " " << ch.roc_first << " " << ch.roc_last << std::endl;
0102         }  // loop over disable channels
0103       }    // loop over the detSetVector
0104     }      // main loop on the map
0105   }        // if file exists
0106 }
0107 
0108 void PixelFEDChannelCollectionMapTestReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0109   edm::ParameterSetDescription desc;
0110   desc.setComment("Reads payloads of type SiPixelFEDChannelContainer");
0111   desc.addUntracked<bool>("printDebug", true);
0112   desc.addUntracked<std::string>("outputFile", "");
0113   descriptions.add("PixelFEDChannelCollectionMapTestReader", desc);
0114 }
0115 
0116 DEFINE_FWK_MODULE(PixelFEDChannelCollectionMapTestReader);