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>
31 #include "TMath.h"
32 
33 class RawPCCProducer : public edm::one::EDProducer<edm::EndLuminosityBlockProducer,
34  edm::one::WatchLuminosityBlocks> {
35  public:
36  explicit RawPCCProducer(const edm::ParameterSet&);
37  ~RawPCCProducer() override;
38 
39  private:
40  void beginLuminosityBlock (edm::LuminosityBlock const& lumiSeg, const edm::EventSetup& iSetup) final;
41  void endLuminosityBlock (edm::LuminosityBlock const& lumiSeg, const edm::EventSetup& iSetup) final;
42  void endLuminosityBlockProduce(edm::LuminosityBlock& lumiSeg, const edm::EventSetup& iSetup) final;
43  void produce (edm::Event& iEvent, const edm::EventSetup& iSetup) final;
44 
45 
47  std::string pccSource_; //input file EDproducer module label
48  std::string prodInst_; //input file product instance
49  std::string takeAverageValue_; //Output average values
50 
51  std::vector<int> modVeto_; //The list of modules to skip in the lumi calc.
52  std::string outputProductName_; //specifies the trigger Rand or ZeroBias
53  std::vector<int> clustersPerBXInput_; //Will fill this with content from PCC
54  std::vector<int> modID_; //vector with Module IDs 1-1 map to bunch x-ing in clusers
55  std::vector<int> events_; //vector with total events at each bxid.
56  std::vector<int> clustersPerBXOutput_; //new vector containing clusters per bxid
57  std::vector<float> corrClustersPerBXOutput_;//new vector containing clusters per bxid with afterglow corrections
58  std::vector<float> errorPerBX_; //Stat error (or number of events)
59  std::vector<int> goodMods_; //The indicies of all the good modules - not vetoed
60  float totalLumi_; //The total raw luminosity from the pixel clusters - not scaled
61  float statErrOnLumi_; //the statistical error on the lumi - large num ie sqrt(N)
62 
65 
66  bool applyCorr_;
67  std::vector<float> correctionScaleFactors_;
68 
69  std::unique_ptr<LumiInfo> outputLumiInfo;
70 
71  std::ofstream csvfile;
72  };
73 
74 //--------------------------------------------------------------------------------------------------
76 {
77  pccSource_ = iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getParameter<std::string>("inputPccLabel");
78  prodInst_ = iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getParameter<std::string>("ProdInst");
79  takeAverageValue_ = iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getUntrackedParameter<std::string>("OutputValue",std::string("Totals"));
80  outputProductName_ = iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getUntrackedParameter<std::string>("outputProductName","alcaLumi");
81  modVeto_ = iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getParameter<std::vector<int>>("modVeto");
82  applyCorr_ = iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getUntrackedParameter<bool>("ApplyCorrections",false);
83  saveCSVFile_ = iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getUntrackedParameter<bool>("saveCSVFile",false);
84  csvOutLabel_ = iConfig.getParameter<edm::ParameterSet>("RawPCCProducerParameters").getUntrackedParameter<std::string>("label",std::string("rawPCC.csv"));
85 
86  edm::InputTag PCCInputTag_(pccSource_, prodInst_);
87  clustersPerBXOutput_.resize(LumiConstants::numBX,0);//vector containing clusters per bxid
88  corrClustersPerBXOutput_.resize(LumiConstants::numBX,0);//vector containing clusters per bxid
89  pccToken=consumes<reco::PixelClusterCounts, edm::InLumi>(PCCInputTag_);
90  produces<LumiInfo, edm::Transition::EndLuminosityBlock>(outputProductName_);
91 
92  if(!applyCorr_){
94  }
95 }
96 
97 //--------------------------------------------------------------------------------------------------
99 }
100 
101 //--------------------------------------------------------------------------------------------------
103 
104 
105 }
106 
107 //--------------------------------------------------------------------------------------------------
109 
110  outputLumiInfo = std::make_unique<LumiInfo>();
111 
112  if(saveCSVFile_){
113  csvfile.open(csvOutLabel_, std::ios_base::app);
114  csvfile<<std::to_string(lumiSeg.run())<<",";
115  csvfile<<std::to_string(lumiSeg.luminosityBlock())<<",";
116  }
117 }
118 
119 //--------------------------------------------------------------------------------------------------
121  totalLumi_=0.0;
122  statErrOnLumi_=0.0;
123 
126 
127  goodMods_.clear();
128 
130  lumiSeg.getByToken(pccToken,pccHandle);
131 
132  const reco::PixelClusterCounts& inputPcc = *(pccHandle.product());
133 
134  modID_ = inputPcc.readModID();
135  events_= inputPcc.readEvents();
136  clustersPerBXInput_ = inputPcc.readCounts();
137 
138  //making list of modules to sum over
139  for (unsigned int i=0;i<modID_.size();i++){
140  if (std::find(modVeto_.begin(),modVeto_.end(), modID_.at(i)) == modVeto_.end()){
141  goodMods_.push_back(i);
142  }
143  }
144 
145  //summing over good modules only
146  for (int bx=0;bx<int(LumiConstants::numBX);bx++){
147  for (unsigned int i=0;i<goodMods_.size();i++){
149 
150  }
151  }
152 
153  if(applyCorr_){
155  iSetup.get<LumiCorrectionsRcd>().get(corrHandle);
156  const LumiCorrections *pccCorrections = corrHandle.product();
157  correctionScaleFactors_ = pccCorrections->getCorrectionsBX();
158  }
159 
160  for (unsigned int i=0;i<clustersPerBXOutput_.size();i++){
161  if (events_.at(i)!=0){
163  }else{
165  }
168  }
169 
170  errorPerBX_.clear();
171  errorPerBX_.assign(events_.begin(),events_.end());
172 
173  if(takeAverageValue_=="Average"){
174  unsigned int NActiveBX=0;
175  for (int bx=0;bx<int(LumiConstants::numBX);bx++){
176  if(events_[bx]>0){
177  NActiveBX++;
178  // Counting where events are will only work
179  // for ZeroBias or AlwaysTrue triggers.
180  // Random triggers will get all BXs.
182  errorPerBX_[bx]=1/TMath::Sqrt(float(events_[bx]));
183  }
184  }
185  if (statErrOnLumi_!=0) {
188  }
189  }
190 
191 
192  outputLumiInfo->setTotalInstLumi(totalLumi_);
193  outputLumiInfo->setTotalInstStatError(statErrOnLumi_);
194 
195  outputLumiInfo->setErrorLumiAllBX(errorPerBX_);
196  outputLumiInfo->setInstLumiAllBX(corrClustersPerBXOutput_);
197 
198  if(saveCSVFile_){
199  csvfile<<std::to_string(totalLumi_);
200 
201  if(totalLumi_>0){
202  for(unsigned int bx=0;bx<LumiConstants::numBX;bx++){
203  csvfile<<","<<std::to_string(corrClustersPerBXOutput_[bx]);
204  }
205  csvfile<<std::endl;
206  } else if (totalLumi_<0) {
207  edm::LogInfo("WARNING")<<"WHY IS LUMI NEGATIVE?!?!?!? "<<totalLumi_;
208  }
209 
210  csvfile.close();
211  }
212 }
213 
214 //--------------------------------------------------------------------------------------------------
217 
218 }
219 
T getParameter(std::string const &) const
std::vector< int > modVeto_
std::vector< float > correctionScaleFactors_
std::unique_ptr< LumiInfo > outputLumiInfo
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
static const unsigned int numBX
Definition: LumiConstants.h:9
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::string prodInst_
std::vector< int > clustersPerBXInput_
LuminosityBlockNumber_t luminosityBlock() const
int iEvent
Definition: GenABIO.cc:230
std::vector< int > const & readCounts() const
std::vector< int > goodMods_
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) final
void put(std::unique_ptr< PROD > product)
Put a new product.
std::vector< int > const & readModID() const
RunNumber_t run() const
std::string outputProductName_
std::string takeAverageValue_
std::vector< int > modID_
std::vector< int > clustersPerBXOutput_
void endLuminosityBlockProduce(edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup) final
void endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &iSetup) final
std::string csvOutLabel_
T const * product() const
Definition: Handle.h:81
const std::vector< float > & getCorrectionsBX() const
~RawPCCProducer() override
std::string pccSource_
T get() const
Definition: EventSetup.h:63
void beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &iSetup) final
edm::EDGetTokenT< reco::PixelClusterCounts > pccToken
RawPCCProducer(const edm::ParameterSet &)
std::vector< float > errorPerBX_
T const * product() const
Definition: ESHandle.h:86
std::ofstream csvfile
std::vector< int > events_
std::vector< int > const & readEvents() const
std::vector< float > corrClustersPerBXOutput_
def move(src, dest)
Definition: eostools.py:510