CMS 3D CMS Logo

HcalLaserEventFilter2012.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <iostream>
4 #include <cmath>
5 #include <sstream>
6 #include <fstream>
7 
8 
9 // user include files
14 
17 #include "zlib.h"
18 
19 using namespace std;
20 #define CHUNK 16384
21 
22 //
23 // constructors and destructor
24 //
26 {
27  verbose_ = ps.getUntrackedParameter<bool>("verbose",false);
28  prefix_ = ps.getUntrackedParameter<std::string>("prefix","");
29  minrun_ = ps.getUntrackedParameter<int>("minrun",-1);
30  maxrun_ = ps.getUntrackedParameter<int>("maxrun",-1);
31  WriteBadToFile_=ps.getUntrackedParameter<bool>("WriteBadToFile",false);
32  if (WriteBadToFile_)
33  outfile_.open("badHcalLaserList_eventfilter.txt");
34  forceFilterTrue_=ps.getUntrackedParameter<bool>("forceFilterTrue",false);
35 
36  minRunInFile=999999; maxRunInFile=1;
37  string eventFileName0=ps.getParameter<string>("eventFileName");
38  string eventFileName = edm::FileInPath(eventFileName0).fullPath();
39 
40  if (verbose_) edm::LogInfo("HcalLaserEventFilter2012") << "HCAL laser event list from file "<<eventFileName;
41  readEventListFile(eventFileName);
42  std::sort(EventList_.begin(), EventList_.end());
43  if (verbose_) edm::LogInfo("HcalLaserEventFilter2012")<<" A total of "<<EventList_.size()<<" listed HCAL laser events found in given run range";
44  if (minrun_==-1 || minrun_<minRunInFile) minrun_=minRunInFile;
45  if (maxrun_==-1 || maxrun_>maxRunInFile) maxrun_=maxRunInFile;
46 }
47 
48 void HcalLaserEventFilter2012::addEventString(const string & eventString)
49 {
50  // Loop through list of bad events, and if run is in allowed range, add bad event to EventList
51  int run=0;
52  unsigned int ls=0;
53  unsigned int event=0;
54  // Check that event list object is in correct form
55  size_t found = eventString.find(":"); // find first colon
56  if (found!=std::string::npos)
57  run=atoi((eventString.substr(0,found)).c_str()); // convert to run
58  else
59  {
60  edm::LogError("HcalLaserEventFilter2012")<<" Unable to parse Event list input '"<<eventString<<"' for run number!";
61  return;
62  }
63  size_t found2 = eventString.find(":",found+1); // find second colon
64  if (found2!=std::string::npos)
65  {
67  ls=atoi((eventString.substr(found+1,(found2-found-1))).c_str()); // convert to ls
68  event=atoi((eventString.substr(found2+1)).c_str()); // convert to event
70  if (ls==0 || event==0) edm::LogWarning("HcalLaserEventFilter2012")<<" Strange lumi, event numbers for input '"<<eventString<<"'";
71  }
72  else
73  {
74  edm::LogError("HcalLaserEventFilter2012")<<"Unable to parse Event list input '"<<eventString<<"' for run number!";
75  return;
76  }
77  // If necessary, check that run is within allowed range
78  if (minrun_>-1 && run<minrun_)
79  {
80  if (verbose_) edm::LogInfo("HcalLaserEventFilter2012") <<"Skipping Event list input '"<<eventString<<"' because it is less than minimum run # "<<minrun_;
81  return;
82  }
83  if (maxrun_>-1 && run>maxrun_)
84  {
85  if (verbose_) edm::LogInfo("HcalLaserEventFilter2012") <<"Skipping Event list input '"<<eventString<<"' because it is greater than maximum run # "<<maxrun_;
86  return;
87  }
88  if (minRunInFile>run) minRunInFile=run;
89  if (maxRunInFile<run) maxRunInFile=run;
90  // Now add event to Event List
91  EventList_.push_back(eventString);
92 }
93 
94 #define LENGTH 0x2000
95 
96 void HcalLaserEventFilter2012::readEventListFile(const string & eventFileName)
97 {
98  gzFile file = gzopen (eventFileName.c_str(), "r");
99  if (! file) {
100  edm::LogError("HcalLaserEventFilter2012")<<" Unable to open event list file "<<eventFileName;
101  return;
102  }
103  string b2;
104  int err;
105  int bytes_read;
106  char buffer[LENGTH];
107  unsigned int i;
108  char * pch;
109 
110  while (true) {
111  bytes_read = gzread (file, buffer, LENGTH - 1);
112  buffer[bytes_read] = '\0';
113  i=0;
114  char* saveptr;
115  pch = strtok_r (buffer,"\n",&saveptr);
116  if (buffer[0] == '\n' ) {
117  addEventString(b2);
118  ++i;
119  } else b2+=pch;
120 
121  while (pch != nullptr)
122  {
123  i+=strlen(pch)+1;
124  if (i>b2.size()) b2= pch;
125  if (i==(LENGTH-1)) {
126  if ((buffer[LENGTH-2] == '\n' )|| (buffer[LENGTH-2] == '\0' )){
127  addEventString(b2);
128  b2="";
129  }
130  } else if (i<LENGTH) {
131  addEventString(b2);
132  }
133  pch = strtok_r (nullptr, "\n", &saveptr);
134  }
135  if (bytes_read < LENGTH - 1) {
136  if (gzeof (file)) break;
137  else {
138  const char * error_string;
139  error_string = gzerror (file, & err);
140  if (err) {
141  edm::LogError("HcalLaserEventFilter2012")<<"Error while reading gzipped file: "<<error_string;
142  return;
143  }
144  }
145  }
146  }
147  gzclose (file);
148  return;
149 }
150 
151 
153 {
154 
155  // do anything here that needs to be done at desctruction time
156  // (e.g. close files, deallocate resources etc.)
157 
158 }
159 
160 
161 //
162 // member functions
163 //
164 
165 // ------------ method called on each new Event ------------
166 bool
168 {
169  int run = iEvent.id().run();
170  // if run is outside filter range, then always return true
171  if (minrun_>-1 && run<minrun_) return true;
172  if (maxrun_>-1 && run>maxrun_) return true;
173 
174  // Okay, now create a string object for this run:ls:event
175  std::stringstream thisevent;
176  thisevent<<run<<":"<<iEvent.luminosityBlock()<<":"<<iEvent.id().event();
177 
178  // Event not found in bad list; it is a good event
179  strVecI it = std::lower_bound(EventList_.begin(), EventList_.end(), thisevent.str());
180  if (it == EventList_.end() || thisevent.str() < *it) return true;
181  // Otherwise, this is a bad event
182  // if verbose, dump out event info
183  // Dump out via cout, or via LogInfo? For now, use cout
184  if (verbose_) std::cout <<prefix_<<thisevent.str()<<std::endl;
185 
186  // To use if we decide on LogInfo:
187  // if (verbose_) edm::LogInfo(prefix_)<<thisevent.str();
188 
189  // Write bad event to file
190  if (WriteBadToFile_)
191  outfile_<<iEvent.id().run()<<":"<<iEvent.luminosityBlock()<<":"<<iEvent.id().event()<<std::endl;
192  if (forceFilterTrue_) return true; // if special input boolean set, always return true, regardless of filter decision
193  return false;
194 }
195 
196 // ------------ method called once each job just after ending the event loop ------------
197 void
199  if (WriteBadToFile_) outfile_.close();
200 }
201 
202 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
203 void
205  //The following says we do not know what parameters are allowed so do no validation
206  // Please change this to state exactly what you do use, even if it is no parameters
208  desc.setUnknown();
209  descriptions.addDefault(desc);
210 }
RunNumber_t run() const
Definition: EventID.h:39
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
HcalLaserEventFilter2012(const edm::ParameterSet &)
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:61
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void addDefault(ParameterSetDescription const &psetDescription)
#define LENGTH
void addEventString(const std::string &eventString)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def ls(path, rec=False)
Definition: eostools.py:349
edm::EventID id() const
Definition: EventBase.h:59
bool filter(edm::Event &, const edm::EventSetup &) override
std::string fullPath() const
Definition: FileInPath.cc:163
void readEventListFile(const std::string &eventFileName)
std::vector< std::string >::iterator strVecI
Definition: event.py:1