CMS 3D CMS Logo

RawPCCProducer.cc
Go to the documentation of this file.
1 
10 #include <string>
11 #include <iostream>
12 #include <fstream>
13 #include <vector>
14 #include <mutex>
15 #include <cmath>
33 
34 class RawPCCProducer : public edm::global::EDProducer<edm::EndLuminosityBlockProducer> {
35 public:
36  explicit RawPCCProducer(const edm::ParameterSet&);
37  ~RawPCCProducer() override;
38 
39 private:
40  void globalEndLuminosityBlockProduce(edm::LuminosityBlock& lumiSeg, const edm::EventSetup& iSetup) const final;
41  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const final;
42 
45  const std::string takeAverageValue_; //Output average values
46 
47  const std::vector<int> modVeto_; //The list of modules to skip in the lumi calc.
48 
51  const bool saveCSVFile_;
52 
53  const bool applyCorr_;
54 };
55 
56 //--------------------------------------------------------------------------------------------------
58  : putToken_{produces<LumiInfo, edm::Transition::EndLuminosityBlock>(
59  iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters")
60  .getUntrackedParameter<std::string>("outputProductName", "alcaLumi"))},
61  takeAverageValue_{iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters")
62  .getUntrackedParameter<std::string>("OutputValue", std::string("Totals"))},
63  modVeto_{
64  iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getParameter<std::vector<int>>("modVeto")},
65  csvOutLabel_{iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters")
66  .getUntrackedParameter<std::string>("label", std::string("rawPCC.csv"))},
67  saveCSVFile_{iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters")
68  .getUntrackedParameter<bool>("saveCSVFile", false)},
69  applyCorr_{iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters")
70  .getUntrackedParameter<bool>("ApplyCorrections", false)} {
71  auto pccSource =
72  iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getParameter<std::string>("inputPccLabel");
73  auto prodInst =
74  iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getParameter<std::string>("ProdInst");
75 
76  pccToken_ = consumes<reco::PixelClusterCounts, edm::InLumi>(edm::InputTag(pccSource, prodInst));
77 }
78 
79 //--------------------------------------------------------------------------------------------------
81 
82 //--------------------------------------------------------------------------------------------------
84 
85 //--------------------------------------------------------------------------------------------------
87  const edm::EventSetup& iSetup) const {
88  float totalLumi = 0.0; //The total raw luminosity from the pixel clusters - not scaled
89  float statErrOnLumi = 0.0; //the statistical error on the lumi - large num ie sqrt(N)
90 
91  //new vector containing clusters per bxid
92  std::vector<int> clustersPerBXOutput(LumiConstants::numBX, 0);
93  //new vector containing clusters per bxid with afterglow corrections
94  std::vector<float> corrClustersPerBXOutput(LumiConstants::numBX, 0);
95 
96  //The indicies of all the good modules - not vetoed
97  std::vector<int> goodMods;
98 
100  lumiSeg.getByToken(pccToken_, pccHandle);
101 
102  const reco::PixelClusterCounts& inputPcc = *(pccHandle.product());
103 
104  //vector with Module IDs 1-1 map to bunch x-ing in clusers
105  auto modID = inputPcc.readModID();
106  //vector with total events at each bxid.
107  auto events = inputPcc.readEvents();
108  auto clustersPerBXInput = inputPcc.readCounts();
109 
110  //making list of modules to sum over
111  for (unsigned int i = 0; i < modID.size(); i++) {
112  if (std::find(modVeto_.begin(), modVeto_.end(), modID.at(i)) == modVeto_.end()) {
113  goodMods.push_back(i);
114  }
115  }
116 
117  //summing over good modules only
118  for (int bx = 0; bx < int(LumiConstants::numBX); bx++) {
119  for (unsigned int i = 0; i < goodMods.size(); i++) {
120  clustersPerBXOutput.at(bx) += clustersPerBXInput.at(goodMods.at(i) * int(LumiConstants::numBX) + bx);
121  }
122  }
123 
124  std::vector<float> correctionScaleFactors;
125  if (applyCorr_) {
127  iSetup.get<LumiCorrectionsRcd>().get(corrHandle);
128  const LumiCorrections* pccCorrections = corrHandle.product();
129  correctionScaleFactors = pccCorrections->getCorrectionsBX();
130  } else {
131  correctionScaleFactors.resize(LumiConstants::numBX, 1.0);
132  }
133 
134  for (unsigned int i = 0; i < clustersPerBXOutput.size(); i++) {
135  if (events.at(i) != 0) {
136  corrClustersPerBXOutput[i] = clustersPerBXOutput[i] * correctionScaleFactors[i];
137  } else {
138  corrClustersPerBXOutput[i] = 0.0;
139  }
140  totalLumi += corrClustersPerBXOutput[i];
141  statErrOnLumi += float(events[i]);
142  }
143 
144  std::vector<float> errorPerBX; //Stat error (or number of events)
145  errorPerBX.assign(events.begin(), events.end());
146 
147  if (takeAverageValue_ == "Average") {
148  unsigned int NActiveBX = 0;
149  for (int bx = 0; bx < int(LumiConstants::numBX); bx++) {
150  if (events[bx] > 0) {
151  NActiveBX++;
152  // Counting where events are will only work
153  // for ZeroBias or AlwaysTrue triggers.
154  // Random triggers will get all BXs.
155  corrClustersPerBXOutput[bx] /= float(events[bx]);
156  errorPerBX[bx] = 1 / sqrt(float(events[bx]));
157  }
158  }
159  if (statErrOnLumi != 0) {
160  totalLumi = totalLumi / statErrOnLumi * float(NActiveBX);
161  statErrOnLumi = 1 / sqrt(statErrOnLumi) * totalLumi;
162  }
163  }
164 
165  LumiInfo outputLumiInfo;
166 
167  outputLumiInfo.setTotalInstLumi(totalLumi);
168  outputLumiInfo.setTotalInstStatError(statErrOnLumi);
169 
170  outputLumiInfo.setErrorLumiAllBX(errorPerBX);
171  outputLumiInfo.setInstLumiAllBX(corrClustersPerBXOutput);
172 
173  if (saveCSVFile_) {
174  std::lock_guard<std::mutex> lock(fileLock_);
175  std::ofstream csfile(csvOutLabel_, std::ios_base::app);
176  csfile << std::to_string(lumiSeg.run()) << ",";
177  csfile << std::to_string(lumiSeg.luminosityBlock()) << ",";
178  csfile << std::to_string(totalLumi);
179 
180  if (totalLumi > 0) {
181  for (unsigned int bx = 0; bx < LumiConstants::numBX; bx++) {
182  csfile << "," << std::to_string(corrClustersPerBXOutput[bx]);
183  }
184  csfile << std::endl;
185  } else if (totalLumi < 0) {
186  edm::LogInfo("WARNING") << "WHY IS LUMI NEGATIVE?!?!?!? " << totalLumi;
187  }
188 
189  csfile.close();
190  }
191  lumiSeg.emplace(putToken_, std::move(outputLumiInfo));
192 }
193 
T getParameter(std::string const &) const
static boost::mutex mutex
Definition: Proxy.cc:11
bool getByToken(EDGetToken token, Handle< PROD > &result) const
const std::string csvOutLabel_
static const unsigned int numBX
Definition: LumiConstants.h:9
void setInstLumiAllBX(std::vector< float > &instLumiByBX)
Definition: LumiInfo.cc:39
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void emplace(EDPutTokenT< PROD > token, Args &&...args)
puts a new product
const std::string takeAverageValue_
LuminosityBlockNumber_t luminosityBlock() const
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< int > const & readCounts() const
T sqrt(T t)
Definition: SSEVec.h:18
void setTotalInstLumi(float totalLumi)
Definition: LumiInfo.h:130
std::vector< int > const & readModID() const
const bool applyCorr_
RunNumber_t run() const
const edm::EDPutTokenT< LumiInfo > putToken_
T const * product() const
Definition: Handle.h:74
const std::vector< float > & getCorrectionsBX() const
std::mutex fileLock_
~RawPCCProducer() override
void setTotalInstStatError(float statError)
Definition: LumiInfo.h:134
const std::vector< int > modVeto_
T get() const
Definition: EventSetup.h:71
void setErrorLumiAllBX(std::vector< float > &errLumiByBX)
Definition: LumiInfo.cc:43
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const final
const bool saveCSVFile_
RawPCCProducer(const edm::ParameterSet &)
T const * product() const
Definition: ESHandle.h:86
std::vector< int > const & readEvents() const
def move(src, dest)
Definition: eostools.py:511
edm::EDGetTokenT< reco::PixelClusterCounts > pccToken_
void globalEndLuminosityBlockProduce(edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup) const final