Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:13

0001 #include "IORawData/CaloPatterns/interface/HcalPatternXMLParser.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "CondFormats/DataRecord/interface/HcalElectronicsMapRcd.h"
0008 #include "CondFormats/HcalObjects/interface/HcalElectronicsMap.h"
0009 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0010 #include <wordexp.h>
0011 #include <cstdio>
0012 
0013 #include <vector>
0014 #include "IORawData/CaloPatterns/interface/HcalFiberPattern.h"
0015 #include "FWCore/Framework/interface/one/EDProducer.h"
0016 
0017 class HcalPatternSource : public edm::one::EDProducer<> {
0018 public:
0019   HcalPatternSource(const edm::ParameterSet& pset);
0020   ~HcalPatternSource() override;
0021   void produce(edm::Event& e, const edm::EventSetup& c) override;
0022 
0023 private:
0024   void loadPatterns(const std::string& patspec);
0025   void loadPatternFile(const std::string& filename);
0026   std::vector<int> bunches_;
0027   std::vector<HcalFiberPattern> patterns_;
0028   int presamples_, samples_;
0029   const edm::ESGetToken<HcalElectronicsMap, HcalElectronicsMapRcd> theHcalElectronicsMapToken_;
0030 };
0031 
0032 HcalPatternSource::HcalPatternSource(const edm::ParameterSet& pset)
0033     : bunches_(pset.getUntrackedParameter<std::vector<int> >("Bunches", std::vector<int>())),
0034       presamples_(pset.getUntrackedParameter<int>("Presamples", 4)),
0035       samples_(pset.getUntrackedParameter<int>("Samples", 10)),
0036       theHcalElectronicsMapToken_(esConsumes()) {
0037   loadPatterns(pset.getUntrackedParameter<std::string>("Patterns"));
0038   produces<HBHEDigiCollection>();
0039   produces<HODigiCollection>();
0040   produces<HFDigiCollection>();
0041 }
0042 
0043 HcalPatternSource::~HcalPatternSource() = default;
0044 
0045 void HcalPatternSource::produce(edm::Event& e, const edm::EventSetup& iSetup) {
0046   if (e.id().event() > bunches_.size())
0047     return;
0048 
0049   // Get HcalElectronicsMap from the Event setup
0050   const auto elecmap = &iSetup.getData(theHcalElectronicsMapToken_);
0051 
0052   auto hbhe = std::make_unique<HBHEDigiCollection>();
0053   auto hf = std::make_unique<HFDigiCollection>();
0054   auto ho = std::make_unique<HODigiCollection>();
0055 
0056   int bc = bunches_[e.id().event() - 1];
0057   for (std::vector<HcalFiberPattern>::iterator i = patterns_.begin(); i != patterns_.end(); i++) {
0058     std::vector<HcalQIESample> samples;
0059     for (int fc = 0; fc < 3; fc++) {
0060       samples = i->getSamples(bc, presamples_, samples_, fc);
0061       HcalElectronicsId eid = i->getId(fc);
0062 
0063       HcalDetId did(elecmap->lookup(eid));
0064 
0065       if (did.null()) {
0066         edm::LogWarning("HcalPatternSource") << "No electronics map match for id " << eid;
0067         continue;
0068       }
0069 
0070       switch (did.subdet()) {
0071         case (HcalBarrel):
0072         case (HcalEndcap):
0073           hbhe->push_back(HBHEDataFrame(did));
0074           hbhe->back().setSize(samples_);
0075           hbhe->back().setPresamples(presamples_);
0076           for (int i = 0; i < samples_; i++)
0077             hbhe->back().setSample(i, samples[i]);
0078           hbhe->back().setReadoutIds(eid);
0079           break;
0080         case (HcalForward):
0081           hf->push_back(HFDataFrame(did));
0082           hf->back().setSize(samples_);
0083           hf->back().setPresamples(presamples_);
0084           for (int i = 0; i < samples_; i++)
0085             hf->back().setSample(i, samples[i]);
0086           hf->back().setReadoutIds(eid);
0087           break;
0088         case (HcalOuter):
0089           ho->push_back(HODataFrame(did));
0090           ho->back().setSize(samples_);
0091           ho->back().setPresamples(presamples_);
0092           for (int i = 0; i < samples_; i++)
0093             ho->back().setSample(i, samples[i]);
0094           ho->back().setReadoutIds(eid);
0095           break;
0096         default:
0097           continue;
0098       }
0099     }
0100   }
0101   hbhe->sort();
0102   ho->sort();
0103   hf->sort();
0104 
0105   e.put(std::move(hbhe));
0106   e.put(std::move(ho));
0107   e.put(std::move(hf));
0108 }
0109 
0110 void HcalPatternSource::loadPatterns(const std::string& patspec) {
0111   wordexp_t p;
0112   char** files;
0113   wordexp(patspec.c_str(), &p, WRDE_NOCMD);  // do not run shell commands!
0114   files = p.we_wordv;
0115   for (unsigned int i = 0; i < p.we_wordc; i++) {
0116     LogDebug("HCAL") << "Reading pattern file '" << files[i] << "'";
0117     loadPatternFile(files[i]);
0118     LogDebug("HCAL") << "Fibers so far " << patterns_.size();
0119   }
0120   wordfree(&p);
0121 }
0122 
0123 void HcalPatternSource::loadPatternFile(const std::string& filename) {
0124   HcalPatternXMLParser parser;
0125   std::string buffer, element;
0126   std::map<std::string, std::string> params;
0127   std::vector<uint32_t> data;
0128   FILE* f = fopen(filename.c_str(), "r");
0129   if (f == nullptr)
0130     return;
0131   else {
0132     char block[4096];
0133     while (!feof(f)) {
0134       int read = fread(block, 1, 4096, f);
0135       buffer.append(block, block + read);
0136     }
0137     fclose(f);
0138   }
0139   if (buffer.find("<?xml") != 0) {
0140     throw cms::Exception("InvalidFormat") << "Not a valid XML file: " << filename;
0141   }
0142   std::string::size_type i = 0, j;
0143   while (buffer.find("<CFGBrick>", i) != std::string::npos) {
0144     i = buffer.find("<CFGBrick>", i);
0145     j = buffer.find("</CFGBrick>", i);
0146     element = "<?xml version='1.0'?>\n";
0147     element.append(buffer, i, j - i);
0148     element.append("</CFGBrick>");
0149     //    LogDebug("HCAL") << element;
0150     params.clear();
0151     data.clear();
0152     parser.parse(element, params, data);
0153     patterns_.push_back(HcalFiberPattern(params, data));
0154     i = j + 5;
0155   }
0156 }
0157 
0158 DEFINE_FWK_MODULE(HcalPatternSource);