Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:42

0001 // system include files
0002 #include <memory>
0003 #include <iostream>
0004 #include <cmath>
0005 #include <sstream>
0006 #include <fstream>
0007 
0008 // user include files
0009 #include "EventFilter/HcalRawToDigi/interface/HcalLaserEventFilter2012.h"
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 #include "FWCore/Framework/interface/Run.h"
0012 #include "FWCore/Framework/interface/LuminosityBlock.h"
0013 
0014 #include "DataFormats/Provenance/interface/EventID.h"
0015 #include "DataFormats/Provenance/interface/RunID.h"
0016 #include "zlib.h"
0017 
0018 using namespace std;
0019 #define CHUNK 16384
0020 
0021 //
0022 // constructors and destructor
0023 //
0024 HcalLaserEventFilter2012::HcalLaserEventFilter2012(const edm::ParameterSet& ps) {
0025   verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
0026   prefix_ = ps.getUntrackedParameter<std::string>("prefix", "");
0027   minrun_ = ps.getUntrackedParameter<int>("minrun", -1);
0028   maxrun_ = ps.getUntrackedParameter<int>("maxrun", -1);
0029   WriteBadToFile_ = ps.getUntrackedParameter<bool>("WriteBadToFile", false);
0030   if (WriteBadToFile_)
0031     outfile_.open("badHcalLaserList_eventfilter.txt");
0032   forceFilterTrue_ = ps.getUntrackedParameter<bool>("forceFilterTrue", false);
0033 
0034   minRunInFile = 999999;
0035   maxRunInFile = 1;
0036   string eventFileName0 = ps.getParameter<string>("eventFileName");
0037   string eventFileName = edm::FileInPath(eventFileName0).fullPath();
0038 
0039   if (verbose_)
0040     edm::LogInfo("HcalLaserEventFilter2012") << "HCAL laser event list from file " << eventFileName;
0041   readEventListFile(eventFileName);
0042   std::sort(EventList_.begin(), EventList_.end());
0043   if (verbose_)
0044     edm::LogInfo("HcalLaserEventFilter2012")
0045         << " A total of " << EventList_.size() << " listed HCAL laser events found in given run range";
0046   if (minrun_ == -1 || minrun_ < minRunInFile)
0047     minrun_ = minRunInFile;
0048   if (maxrun_ == -1 || maxrun_ > maxRunInFile)
0049     maxrun_ = maxRunInFile;
0050 }
0051 
0052 void HcalLaserEventFilter2012::addEventString(const string& eventString) {
0053   // Loop through list of bad events, and if run is in allowed range, add bad event to EventList
0054   int run = 0;
0055   unsigned int ls = 0;
0056   unsigned int event = 0;
0057   // Check that event list object is in correct form
0058   size_t found = eventString.find(':');  // find first colon
0059   if (found != std::string::npos)
0060     run = atoi((eventString.substr(0, found)).c_str());  // convert to run
0061   else {
0062     edm::LogError("HcalLaserEventFilter2012")
0063         << "  Unable to parse Event list input '" << eventString << "' for run number!";
0064     return;
0065   }
0066   size_t found2 = eventString.find(':', found + 1);  // find second colon
0067   if (found2 != std::string::npos) {
0068     /// Some event numbers are less than 0?  \JetHT\Run2012C-v1\RAW:201278:2145:-2130281065  -- due to events being dumped out as ints, not uints!
0069     ls = atoi((eventString.substr(found + 1, (found2 - found - 1))).c_str());  // convert to ls
0070     event = atoi((eventString.substr(found2 + 1)).c_str());                    // convert to event
0071     /// Some event numbers are less than 0?  \JetHT\Run2012C-v1\RAW:201278:2145:-2130281065
0072     if (ls == 0 || event == 0)
0073       edm::LogWarning("HcalLaserEventFilter2012") << "  Strange lumi, event numbers for input '" << eventString << "'";
0074   } else {
0075     edm::LogError("HcalLaserEventFilter2012")
0076         << "Unable to parse Event list input '" << eventString << "' for run number!";
0077     return;
0078   }
0079   // If necessary, check that run is within allowed range
0080   if (minrun_ > -1 && run < minrun_) {
0081     if (verbose_)
0082       edm::LogInfo("HcalLaserEventFilter2012")
0083           << "Skipping Event list input '" << eventString << "' because it is less than minimum run # " << minrun_;
0084     return;
0085   }
0086   if (maxrun_ > -1 && run > maxrun_) {
0087     if (verbose_)
0088       edm::LogInfo("HcalLaserEventFilter2012")
0089           << "Skipping Event list input '" << eventString << "' because it is greater than maximum run # " << maxrun_;
0090     return;
0091   }
0092   if (minRunInFile > run)
0093     minRunInFile = run;
0094   if (maxRunInFile < run)
0095     maxRunInFile = run;
0096   // Now add event to Event List
0097   EventList_.push_back(eventString);
0098 }
0099 
0100 #define LENGTH 0x2000
0101 
0102 void HcalLaserEventFilter2012::readEventListFile(const string& eventFileName) {
0103   gzFile file = gzopen(eventFileName.c_str(), "r");
0104   if (!file) {
0105     edm::LogError("HcalLaserEventFilter2012") << "  Unable to open event list file " << eventFileName;
0106     return;
0107   }
0108   string b2;
0109   int err;
0110   int bytes_read;
0111   char buffer[LENGTH];
0112   unsigned int i;
0113   char* pch;
0114 
0115   while (true) {
0116     bytes_read = gzread(file, buffer, LENGTH - 1);
0117     buffer[bytes_read] = '\0';
0118     i = 0;
0119     char* saveptr;
0120     pch = strtok_r(buffer, "\n", &saveptr);
0121     if (buffer[0] == '\n') {
0122       addEventString(b2);
0123       ++i;
0124     } else
0125       b2 += pch;
0126 
0127     while (pch != nullptr) {
0128       i += strlen(pch) + 1;
0129       if (i > b2.size())
0130         b2 = pch;
0131       if (i == (LENGTH - 1)) {
0132         if ((buffer[LENGTH - 2] == '\n') || (buffer[LENGTH - 2] == '\0')) {
0133           addEventString(b2);
0134           b2 = "";
0135         }
0136       } else if (i < LENGTH) {
0137         addEventString(b2);
0138       }
0139       pch = strtok_r(nullptr, "\n", &saveptr);
0140     }
0141     if (bytes_read < LENGTH - 1) {
0142       if (gzeof(file))
0143         break;
0144       else {
0145         const char* error_string;
0146         error_string = gzerror(file, &err);
0147         if (err) {
0148           edm::LogError("HcalLaserEventFilter2012") << "Error while reading gzipped file:  " << error_string;
0149           return;
0150         }
0151       }
0152     }
0153   }
0154   gzclose(file);
0155   return;
0156 }
0157 
0158 HcalLaserEventFilter2012::~HcalLaserEventFilter2012() {
0159   // do anything here that needs to be done at desctruction time
0160   // (e.g. close files, deallocate resources etc.)
0161 }
0162 
0163 //
0164 // member functions
0165 //
0166 
0167 // ------------ method called on each new Event  ------------
0168 bool HcalLaserEventFilter2012::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0169   int run = iEvent.id().run();
0170   // if run is outside filter range, then always return true
0171   if (minrun_ > -1 && run < minrun_)
0172     return true;
0173   if (maxrun_ > -1 && run > maxrun_)
0174     return true;
0175 
0176   // Okay, now create a string object for this run:ls:event
0177   std::stringstream thisevent;
0178   thisevent << run << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event();
0179 
0180   // Event not found in bad list; it is a good event
0181   strVecI it = std::lower_bound(EventList_.begin(), EventList_.end(), thisevent.str());
0182   if (it == EventList_.end() || thisevent.str() < *it)
0183     return true;
0184   // Otherwise, this is a bad event
0185   // if verbose, dump out event info
0186   // Dump out via cout, or via LogInfo?  For now, use cout
0187   if (verbose_)
0188     std::cout << prefix_ << thisevent.str() << std::endl;
0189 
0190   // To use if we decide on LogInfo:
0191   // if (verbose_) edm::LogInfo(prefix_)<<thisevent.str();
0192 
0193   // Write bad event to file
0194   if (WriteBadToFile_)
0195     outfile_ << iEvent.id().run() << ":" << iEvent.luminosityBlock() << ":" << iEvent.id().event() << std::endl;
0196   if (forceFilterTrue_)
0197     return true;  // if special input boolean set, always return true, regardless of filter decision
0198   return false;
0199 }
0200 
0201 // ------------ method called once each job just after ending the event loop  ------------
0202 void HcalLaserEventFilter2012::endJob() {
0203   if (WriteBadToFile_)
0204     outfile_.close();
0205 }
0206 
0207 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0208 void HcalLaserEventFilter2012::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0209   //The following says we do not know what parameters are allowed so do no validation
0210   // Please change this to state exactly what you do use, even if it is no parameters
0211   edm::ParameterSetDescription desc;
0212   desc.setUnknown();
0213   descriptions.addDefault(desc);
0214 }
0215 DEFINE_FWK_MODULE(HcalLaserEventFilter2012);