CMS 3D CMS Logo

LumiProducerFromBrilcalc.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoLuminosity/LumiProducer
4 // Class: LumiProducerFromBrilcalc
5 //
13 //
14 // Original Author: Paul Lujan
15 // Created: Fri, 28 Feb 2020 16:32:38 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <sstream>
22 #include <fstream>
23 
24 // user include files
27 
30 
33 
35 
37 
38 //
39 // class declaration
40 //
41 
43 public:
45  ~LumiProducerFromBrilcalc() override = default;
46 
47  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
48 
49 private:
50  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
51 
52  // ----------member data ---------------------------
54  const bool throwIfNotFound_;
55  const bool doBunchByBunch_;
56  std::map<std::pair<int, int>, std::pair<float, float> > lumiData_;
57 };
58 
59 //
60 // constructor
61 //
62 
64  : lumiFile_(iConfig.getParameter<std::string>("lumiFile")),
65  throwIfNotFound_(iConfig.getParameter<bool>("throwIfNotFound")),
66  doBunchByBunch_(iConfig.getParameter<bool>("doBunchByBunch")) {
67  //register your products
68  produces<LumiInfo>("brilcalc");
69 
70  //now do what ever other initialization is needed
71  if (doBunchByBunch_) {
72  throw cms::Exception("LumiProducerFromBrilcalc")
73  << "Sorry, bunch-by-bunch luminosity is not yet supported. Please bug your friendly lumi expert!";
74  }
75 
76  // Read luminosity data and store it in lumiData_.
77  edm::LogInfo("LumiProducerFromBrilcalc") << "Reading luminosity data from " << lumiFile_ << "...one moment...";
78  std::ifstream lumiFile(lumiFile_);
79  if (!lumiFile.is_open()) {
80  throw cms::Exception("LumiProducerFromBrilcalc") << "Failed to open input luminosity file " << lumiFile_;
81  }
82 
83  int nLS = 0;
85  while (true) {
86  std::getline(lumiFile, line);
87  if (lumiFile.eof() || lumiFile.fail())
88  break;
89  if (line.empty())
90  continue; // skip blank lines
91  if (line.at(0) == '#')
92  continue; // skip comment lines
93 
94  // Break into fields. These should be, in order: run:fill, brills:cmsls, time, beam status, beam energy,
95  // delivered lumi, recorded lumi, average pileup, source.
96  std::stringstream ss(line);
97  std::string field;
98  std::vector<std::string> fields;
99 
100  while (std::getline(ss, field, ','))
101  fields.push_back(field);
102 
103  if (fields.size() != 9) {
104  edm::LogWarning("LumiProducerFromBrilcalc") << "Malformed line in csv file: " << line;
105  continue;
106  }
107 
108  // Extract the run number from fields[0] and the lumisection number from fields[1]. Fortunately since the
109  // thing we want is before the colon, we don't have to split them again.
110  int run, ls;
111  std::stringstream runfill(fields[0]);
112  runfill >> run;
113  std::stringstream lsls(fields[1]);
114  lsls >> ls;
115 
116  // Extract the delivered and recorded lumi from fields[5] and fields[6].
117  float lumiDel, lumiRec, dtFrac;
118  std::stringstream lumiDelString(fields[5]);
119  lumiDelString >> lumiDel;
120  std::stringstream lumiRecString(fields[6]);
121  lumiRecString >> lumiRec;
122 
123  // Calculate the deadtime fraction
124  dtFrac = 1.0 - lumiRec / lumiDel;
125 
126  // Finally, store it all.
127  lumiData_[std::make_pair(run, ls)] = std::make_pair(lumiDel, dtFrac);
128  nLS++;
129  }
130  edm::LogInfo("LumiProducerFromBrilcalc") << "Read " << nLS << " lumisections from " << lumiFile_;
131  lumiFile.close();
132 }
133 
134 //
135 // member functions
136 //
137 
138 // ------------ method called to produce the data ------------
141  const edm::EventSetup& iSetup) const {
142  std::vector<float> bxlumi(3564, 0);
143  std::pair<int, int> runls = std::make_pair(iEvent.run(), iEvent.luminosityBlock());
144  if (lumiData_.count(runls) == 1) {
145  // if we have data for this run/LS, put it in the event
146  LogDebug("LumiProducerFromBrilcalc") << "Filling for run " << runls.first << " ls " << runls.second
147  << " with delivered " << lumiData_.at(runls).first << " dt "
148  << lumiData_.at(runls).second;
149  iEvent.put(std::make_unique<LumiInfo>(lumiData_.at(runls).second, bxlumi, lumiData_.at(runls).first), "brilcalc");
150  } else {
151  if (throwIfNotFound_) {
152  throw cms::Exception("LumiProducerFromBrilcalc")
153  << "Failed to find luminosity for run " << runls.first << " LS " << runls.second;
154  } else {
155  // just put in zeroes
156  edm::LogWarning("LumiProducerFromBrilcalc")
157  << "Failed to find luminosity for run " << runls.first << " ls " << runls.second;
158  iEvent.put(std::make_unique<LumiInfo>(0, bxlumi, 0), "brilcalc");
159  }
160  }
161 }
162 
163 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
165  // Allowed parameters
167  desc.add<std::string>("lumiFile");
168  desc.add<bool>("throwIfNotFound", false);
169  desc.add<bool>("doBunchByBunch", false);
170  descriptions.addDefault(desc);
171 }
172 
173 //define this as a plug-in
~LumiProducerFromBrilcalc() override=default
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
LumiProducerFromBrilcalc(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::map< std::pair< int, int >, std::pair< float, float > > lumiData_
Log< level::Info, false > LogInfo
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Log< level::Warning, false > LogWarning
#define LogDebug(id)