CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
CorrPCCProducer Class Reference
Inheritance diagram for CorrPCCProducer:
edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::EndLuminosityBlockProducer, edm::one::WatchLuminosityBlocks > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 CorrPCCProducer (const edm::ParameterSet &)
 
 ~CorrPCCProducer () override
 
- Public Member Functions inherited from edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::EndLuminosityBlockProducer, edm::one::WatchLuminosityBlocks >
 EDProducer ()=default
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
virtual ~ProducerBase () noexcept(false)
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void beginLuminosityBlock (edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &iSetup) final
 
void beginRun (edm::Run const &runSeg, const edm::EventSetup &iSetup) final
 
void calculateCorrections (std::vector< float >, std::vector< float > &, float &)
 
void endJob () final
 
void endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &iSetup) final
 
void endLuminosityBlockProduce (edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup) final
 
void endRun (edm::Run const &runSeg, const edm::EventSetup &iSetup) final
 
void endRunProduce (edm::Run &runSeg, const edm::EventSetup &iSetup) final
 
void estimateType1Frac (std::vector< float >, float &)
 
void evaluateCorrectionResiduals (std::vector< float >)
 
float getMaximum (std::vector< float >)
 
void makeCorrectionTemplate ()
 
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) final
 
void resetBlock ()
 

Private Attributes

unsigned int approxLumiBlockSize_
 
std::vector< float > correctionScaleFactors_
 
std::vector< float > correctionTemplate_
 
TH1F * corrlumiAvg_h
 
unsigned int countLumi_
 
std::vector< float > errOnLumiByBX_
 
std::vector< float > events_
 
TFile * histoFile
 
TList * hlist
 
unsigned int iBlock =0
 
TH1F * lumiAvg_h
 
std::map< std::pair< unsigned int, unsigned int >, unsigned int > lumiInfoCounter
 
std::map< std::pair< unsigned int, unsigned int >, LumiInfo * > lumiInfoMap
 
std::map< std::pair< unsigned int, unsigned int >, LumiInfo * >::iterator lumiInfoMapIterator
 
std::map< unsigned int, LumiInfo * > lumiInfoMapPerLS
 
edm::EDGetTokenT< LumiInfolumiInfoToken
 
std::vector< unsigned int > lumiSections
 
float mean_type1_residual
 
float mean_type1_residual_unc
 
float mean_type2_residual
 
float mean_type2_residual_unc
 
unsigned int minimumNumberOfEvents
 
unsigned int nTrain
 
float overallCorrection_
 
LumiCorrectionspccCorrections
 
std::string pccSrc_
 
edm::Service< cond::service::PoolDBOutputServicepoolDbService
 
std::string prodInst_
 
std::vector< float > rawlumiBX_
 
TH1F * scaleFactorAvg_h
 
unsigned int thisLS
 
std::vector< float > totalLumiByBX_
 
std::vector< float > totalLumiByBXAvg_
 
float type1Frac
 
TGraphErrors * type1FracGraph
 
TH1F * type1FracHist
 
TGraphErrors * type1resGraph
 
TH1F * type1resHist
 
double type2_a_
 
double type2_b_
 
TGraphErrors * type2resGraph
 
TH1F * type2resHist
 

Additional Inherited Members

- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description


class: CorrPCCProducer.cc

description: Computes the type 1 and type 2 corrections to the luminosity type 1 - first (spillover from previous BXs real clusters) type 2 - after (comes from real activation)

authors:Sam Higginbotham (shigg.nosp@m.inb@.nosp@m.cern..nosp@m.ch) and Chris Palmer (capal.nosp@m.mer@.nosp@m.cern..nosp@m.ch)


Definition at line 48 of file CorrPCCProducer.cc.

Constructor & Destructor Documentation

CorrPCCProducer::CorrPCCProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 124 of file CorrPCCProducer.cc.

References approxLumiBlockSize_, correctionScaleFactors_, correctionTemplate_, countLumi_, events_, edm::ParameterSet::getParameter(), histoFile, lumiInfoToken, makeCorrectionTemplate(), minimumNumberOfEvents, LumiConstants::numBX, pccSrc_, prodInst_, resetBlock(), totalLumiByBX_, totalLumiByBXAvg_, type1FracGraph, type1resGraph, type2_a_, type2_b_, and type2resGraph.

125 {
126  pccSrc_ = iConfig.getParameter<edm::ParameterSet>("CorrPCCProducerParameters").getParameter<std::string>("inLumiObLabel");
127  prodInst_ = iConfig.getParameter<edm::ParameterSet>("CorrPCCProducerParameters").getParameter<std::string>("ProdInst");
128  approxLumiBlockSize_=iConfig.getParameter<edm::ParameterSet>("CorrPCCProducerParameters").getParameter<int>("approxLumiBlockSize");
129  type2_a_ = iConfig.getParameter<edm::ParameterSet>("CorrPCCProducerParameters").getParameter<double>("type2_a");
130  type2_b_ = iConfig.getParameter<edm::ParameterSet>("CorrPCCProducerParameters").getParameter<double>("type2_b");
131  countLumi_=0;
133 
139 
140  resetBlock();
141 
143 
144  edm::InputTag inputPCCTag_(pccSrc_, prodInst_);
145 
146  lumiInfoToken=consumes<LumiInfo, edm::InLumi>(inputPCCTag_);
147 
148  histoFile = new TFile("CorrectionHisto.root","RECREATE");
149 
150  type1FracGraph= new TGraphErrors();
151  type1resGraph = new TGraphErrors();
152  type2resGraph = new TGraphErrors();
153  type1FracGraph->SetName("Type1Fraction");
154  type1resGraph->SetName("Type1Res");
155  type2resGraph->SetName("Type2Res");
156  type1FracGraph->GetYaxis()->SetTitle("Type 1 Fraction");
157  type1resGraph->GetYaxis()->SetTitle("Type 1 Residual");
158  type2resGraph->GetYaxis()->SetTitle("Type 2 Residual");
159  type1FracGraph->GetXaxis()->SetTitle("Unique LS ID");
160  type1resGraph->GetXaxis()->SetTitle("Unique LS ID");
161  type2resGraph->GetXaxis()->SetTitle("Unique LS ID");
162  type1FracGraph->SetMarkerStyle(8);
163  type1resGraph->SetMarkerStyle(8);
164  type2resGraph->SetMarkerStyle(8);
165 }
T getParameter(std::string const &) const
std::vector< float > events_
edm::EDGetTokenT< LumiInfo > lumiInfoToken
static const unsigned int numBX
Definition: LumiConstants.h:9
unsigned int countLumi_
std::string pccSrc_
unsigned int minimumNumberOfEvents
std::vector< float > totalLumiByBXAvg_
TGraphErrors * type1resGraph
std::string prodInst_
std::vector< float > correctionScaleFactors_
std::vector< float > correctionTemplate_
TGraphErrors * type2resGraph
std::vector< float > totalLumiByBX_
TGraphErrors * type1FracGraph
unsigned int approxLumiBlockSize_
void makeCorrectionTemplate()
CorrPCCProducer::~CorrPCCProducer ( )
override

Definition at line 168 of file CorrPCCProducer.cc.

168  {
169 }

Member Function Documentation

void CorrPCCProducer::beginLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
const edm::EventSetup iSetup 
)
finalprivate

Definition at line 350 of file CorrPCCProducer.cc.

References countLumi_.

350  {
351  countLumi_++;
352 }
unsigned int countLumi_
void CorrPCCProducer::beginRun ( edm::Run const &  runSeg,
const edm::EventSetup iSetup 
)
finalprivate

Definition at line 347 of file CorrPCCProducer.cc.

347  {
348 }
void CorrPCCProducer::calculateCorrections ( std::vector< float >  uncorrected,
std::vector< float > &  correctionScaleFactors_,
float &  overallCorrection_ 
)
private

Definition at line 284 of file CorrPCCProducer.cc.

References correctionTemplate_, estimateType1Frac(), evaluateCorrectionResiduals(), getMaximum(), mps_fire::i, SiStripPI::max, LumiConstants::numBX, electronIdCutBased_cfi::threshold, and type1Frac.

Referenced by endRunProduce().

284  {
285  type1Frac = 0;
286 
287  int nTrials=4;
288 
289  for(int trial=0;trial<nTrials;trial++){
290  estimateType1Frac(uncorrected, type1Frac);
291  edm::LogInfo("INFO") <<"type 1 fraction after iteration "<<trial<<" is "<<type1Frac;
292  }
293 
294  //correction should never be negative
295  type1Frac=std::max(0.0,(double)type1Frac);
296 
297  std::vector<float> corrected_tmp_;
298  for(size_t i=0; i<uncorrected.size(); i++){
299  corrected_tmp_.push_back(uncorrected.at(i));
300  }
301 
302 
303  //Apply all corrections
304  for(size_t i=0; i<LumiConstants::numBX-1; i++){
305  // type 1 - first (spillover from previous BXs real clusters)
306  float bin_i = corrected_tmp_.at(i);
307  corrected_tmp_.at(i+1)=corrected_tmp_.at(i+1)-type1Frac*bin_i;
308 
309  // type 2 - after (comes from real activation)
310  bin_i = corrected_tmp_.at(i);
311  for(size_t j=i+1; j<i+LumiConstants::numBX-1; j++){
312  if(j<LumiConstants::numBX){
313  corrected_tmp_.at(j)=corrected_tmp_.at(j)-bin_i*correctionTemplate_.at(j-i);
314  }else{
315  corrected_tmp_.at(j-LumiConstants::numBX) = corrected_tmp_.at(j-LumiConstants::numBX)-bin_i*correctionTemplate_.at(j-i);
316  }
317  }
318  }
319 
320  evaluateCorrectionResiduals(corrected_tmp_);
321 
322  float integral_uncorr_clusters=0;
323  float integral_corr_clusters = 0;
324 
325  //Calculate Per-BX correction factor and overall correction factor
326  float lumiMax = getMaximum(corrected_tmp_);
327  float threshold = lumiMax*0.2;
328  for (size_t ibx=0; ibx<corrected_tmp_.size(); ibx++){
329  if(corrected_tmp_.at(ibx)>threshold){
330  integral_uncorr_clusters+=uncorrected.at(ibx);
331  integral_corr_clusters+=corrected_tmp_.at(ibx);
332  }
333  if(corrected_tmp_.at(ibx)!=0.0&&uncorrected.at(ibx)!=0.0){
334  correctionScaleFactors_.at(ibx) = corrected_tmp_.at(ibx)/uncorrected.at(ibx);
335  }else{
336  correctionScaleFactors_.at(ibx) = 0.0;
337  }
338  }
339 
340  overallCorrection_ = integral_corr_clusters/integral_uncorr_clusters;
341 }
static const unsigned int numBX
Definition: LumiConstants.h:9
void evaluateCorrectionResiduals(std::vector< float >)
std::vector< float > correctionScaleFactors_
std::vector< float > correctionTemplate_
float getMaximum(std::vector< float >)
void estimateType1Frac(std::vector< float >, float &)
void CorrPCCProducer::endJob ( void  )
finalprivatevirtual

Reimplemented from edm::one::EDProducerBase.

Definition at line 620 of file CorrPCCProducer.cc.

References DEFINE_FWK_MODULE, histoFile, type1FracGraph, type1resGraph, and type2resGraph.

Referenced by o2olib.O2ORunMgr::executeJob().

620  {
621  histoFile->cd();
622  type1FracGraph->Write();
623  type1resGraph->Write();
624  type2resGraph->Write();
625  histoFile->Write();
626  histoFile->Close();
627  delete type1FracGraph;
628  delete type1resGraph;
629  delete type2resGraph;
630 }
TGraphErrors * type1resGraph
TGraphErrors * type2resGraph
TGraphErrors * type1FracGraph
void CorrPCCProducer::endLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
const edm::EventSetup iSetup 
)
finalprivate

Definition at line 356 of file CorrPCCProducer.cc.

References errOnLumiByBX_, events_, edm::LuminosityBlock::getByToken(), LumiInfo::getErrorLumiAllBX(), LumiInfo::getInstLumiAllBX(), lumiInfoMapPerLS, lumiInfoToken, edm::LuminosityBlockBase::luminosityBlock(), lumiSections, minimumNumberOfEvents, LumiConstants::numBX, edm::Handle< T >::product(), thisLS, and totalLumiByBX_.

356  {
357 
358  thisLS = lumiSeg.luminosityBlock();
359 
360  edm::Handle<LumiInfo> PCCHandle;
361  lumiSeg.getByToken(lumiInfoToken,PCCHandle);
362 
363  const LumiInfo& inLumiOb = *(PCCHandle.product());
364 
365  errOnLumiByBX_= inLumiOb.getErrorLumiAllBX();
366 
367  unsigned int totalEvents=0;
368  for(unsigned int bx=0;bx<LumiConstants::numBX;bx++){
369  totalEvents+=errOnLumiByBX_[bx];
370  }
371 
372  if(totalEvents<minimumNumberOfEvents){
373  edm::LogInfo("INFO") <<"number of events in this LS is too few "<<totalEvents;
374  //return;
375  }
376  else{
377 
378  edm::LogInfo("INFO") <<"Skipping Lumisection "<<thisLS;
379  }
380 
382  totalLumiByBX_=inLumiOb.getInstLumiAllBX();
383  events_=inLumiOb.getErrorLumiAllBX();
384  lumiInfoMapPerLS[thisLS]->setInstLumiAllBX(totalLumiByBX_);
385  lumiInfoMapPerLS[thisLS]->setErrorLumiAllBX(events_);
386  lumiSections.push_back(thisLS);
387 }
std::vector< float > events_
edm::EDGetTokenT< LumiInfo > lumiInfoToken
static const unsigned int numBX
Definition: LumiConstants.h:9
const std::vector< float > & getErrorLumiAllBX() const
Definition: LumiInfo.h:104
unsigned int minimumNumberOfEvents
const std::vector< float > & getInstLumiAllBX() const
Definition: LumiInfo.h:100
std::vector< unsigned int > lumiSections
T const * product() const
Definition: Handle.h:81
std::vector< float > totalLumiByBX_
std::vector< float > errOnLumiByBX_
unsigned int thisLS
std::map< unsigned int, LumiInfo * > lumiInfoMapPerLS
void CorrPCCProducer::endLuminosityBlockProduce ( edm::LuminosityBlock lumiSeg,
const edm::EventSetup iSetup 
)
finalprivate

Definition at line 394 of file CorrPCCProducer.cc.

394  {
395 }
void CorrPCCProducer::endRun ( edm::Run const &  runSeg,
const edm::EventSetup iSetup 
)
finalprivate

Definition at line 390 of file CorrPCCProducer.cc.

390  {
391 }
void CorrPCCProducer::endRunProduce ( edm::Run runSeg,
const edm::EventSetup iSetup 
)
finalprivate

Definition at line 398 of file CorrPCCProducer.cc.

References approxLumiBlockSize_, calculateCorrections(), correctionScaleFactors_, corrlumiAvg_h, errOnLumiByBX_, events_, plotBeamSpotDB::first, objects.autophobj::float, histoFile, mps_fire::i, iBlock, edm::RunBase::id(), edm::Service< T >::isAvailable(), lumiAvg_h, lumiInfoCounter, lumiInfoMap, lumiInfoMapIterator, lumiInfoMapPerLS, lumiSections, mean_type1_residual, mean_type1_residual_unc, mean_type2_residual, mean_type2_residual_unc, LumiConstants::numBX, overallCorrection_, pccCorrections, poolDbService, rawlumiBX_, resetBlock(), edm::RunID::run(), edm::RunBase::run(), scaleFactorAvg_h, LumiCorrections::setCorrectionsBX(), LumiCorrections::setOverallCorrection(), LumiCorrections::setType1Fraction(), LumiCorrections::setType1Residual(), LumiCorrections::setType2Residual(), thisLS, totalLumiByBX_, totalLumiByBXAvg_, type1Frac, type1FracGraph, type1FracHist, type1resGraph, type1resHist, type2resGraph, type2resHist, and cond::service::PoolDBOutputService::writeOne().

398  {
399  if(lumiSections.empty()){
400  return;
401  }
402 
403  std::sort (lumiSections.begin(), lumiSections.end());
404 
405  edm::LogInfo("INFO") <<"Number of Lumisections "<<lumiSections.size()<<" in run "<<runSeg.run() ;
406 
407  //Determining integer number of blocks
408  float nBlocks_f = float(lumiSections.size())/approxLumiBlockSize_;
409  unsigned int nBlocks=1;
410  if(nBlocks_f>1){
411  if(nBlocks_f - lumiSections.size()/approxLumiBlockSize_ < 0.5){
412  nBlocks=lumiSections.size()/approxLumiBlockSize_;
413  } else {
414  nBlocks=lumiSections.size()/approxLumiBlockSize_+1;
415  }
416  }
417 
418  float nLSPerBlock = float(lumiSections.size())/nBlocks;
419 
420  std::vector<std::pair<unsigned int, unsigned int>> lsKeys;
421  lsKeys.clear();
422 
423  //Constructing nBlocks IOVs
424  for(unsigned iKey=0; iKey<nBlocks; iKey++){
425  lsKeys.push_back(std::make_pair(lumiSections[(unsigned int)(iKey*nLSPerBlock)], lumiSections[(unsigned int)((iKey+1)*nLSPerBlock)-1]));
426  }
427 
428  lsKeys[0].first=1;
429 
430  for(unsigned int lumiSection=0; lumiSection<lumiSections.size(); lumiSection++){
431  thisLS=lumiSections[lumiSection];
432 
433  std::pair<unsigned int,unsigned int> lsKey;
434 
435  bool foundKey=false;
436 
437  for(unsigned iKey=0; iKey<nBlocks; iKey++){
438  if( (thisLS >= lsKeys[iKey].first) && (thisLS<= lsKeys[iKey].second) ){
439  lsKey=lsKeys[iKey];
440  foundKey=true;
441  break;
442  }
443  }
444 
445  if(!foundKey){
446  edm::LogInfo("WARNING")<<"Didn't find key "<<thisLS;
447  continue;
448  }
449 
450  if(lumiInfoMap.count( lsKey ) == 0 ) {
451  lumiInfoMap[lsKey] = new LumiInfo();
452  }
453 
454  //Sum all lumi in IOV of lsKey
455  totalLumiByBX_=lumiInfoMap[lsKey]->getInstLumiAllBX();
456  events_ =lumiInfoMap[lsKey]->getErrorLumiAllBX();
457 
458  rawlumiBX_ =lumiInfoMapPerLS[thisLS]->getInstLumiAllBX();
459  errOnLumiByBX_=lumiInfoMapPerLS[thisLS]->getErrorLumiAllBX();
460 
461  for(unsigned int bx=0;bx<LumiConstants::numBX;bx++){
462  totalLumiByBX_[bx]+=rawlumiBX_[bx];
463  events_[bx]+=errOnLumiByBX_[bx];
464  }
465  lumiInfoMap[lsKey]->setInstLumiAllBX(totalLumiByBX_);
466  lumiInfoMap[lsKey]->setErrorLumiAllBX(events_);
467  lumiInfoCounter[lsKey]++;
468  }
469 
470 
471  cond::Time_t thisIOV = 1;
472 
473  char *histname1 = new char[100];
474  char *histname2 = new char[100];
475  char *histname3 = new char[100];
476  char *histTitle1 = new char[100];
477  char *histTitle2 = new char[100];
478  char *histTitle3 = new char[100];
479  sprintf(histTitle1,"Type1Fraction_%d",runSeg.run());
480  sprintf(histTitle2,"Type1Res_%d",runSeg.run());
481  sprintf(histTitle3,"Type2Res_%d",runSeg.run());
482  type1FracHist = new TH1F(histTitle1,histTitle1,1000,-0.5,0.5);
483  type1resHist = new TH1F(histTitle2,histTitle2,4000,-0.2,0.2);
484  type2resHist = new TH1F(histTitle3,histTitle3,4000,-0.2,0.2);
485  delete[] histTitle1;
486  delete[] histTitle2;
487  delete[] histTitle3;
488 
490 
491  totalLumiByBX_=lumiInfoMapIterator->second->getInstLumiAllBX();
492  events_=lumiInfoMapIterator->second->getErrorLumiAllBX();
493 
494  if(events_.empty()){
495  continue;
496  }
497 
499  thisIOV = (cond::Time_t)(lu.value());
500 
501  sprintf(histname1, "CorrectedLumiAvg_%d_%d_%d_%d",runSeg.run(),iBlock,lumiInfoMapIterator->first.first,lumiInfoMapIterator->first.second);
502  sprintf(histname2, "ScaleFactorsAvg_%d_%d_%d_%d",runSeg.run(),iBlock,lumiInfoMapIterator->first.first,lumiInfoMapIterator->first.second);
503  sprintf(histname3, "RawLumiAvg_%d_%d_%d_%d",runSeg.run(),iBlock,lumiInfoMapIterator->first.first,lumiInfoMapIterator->first.second);
504 
505  corrlumiAvg_h = new TH1F(histname1,"",LumiConstants::numBX,1,LumiConstants::numBX);
506  scaleFactorAvg_h= new TH1F(histname2,"",LumiConstants::numBX,1,LumiConstants::numBX);
507  lumiAvg_h = new TH1F(histname3,"",LumiConstants::numBX,1,LumiConstants::numBX);
508 
509 
510  //Averaging by the number of events
511  for(unsigned int i=0;i<LumiConstants::numBX;i++){
512  if(events_.at(i)!=0){
514  }else{
515  totalLumiByBXAvg_[i]=0.0;
516  }
517  }
518 
520 
521  for(unsigned int bx=0;bx<LumiConstants::numBX;bx++){
522  corrlumiAvg_h->SetBinContent(bx,totalLumiByBXAvg_[bx]*correctionScaleFactors_[bx]);
523  if(events_.at(bx)!=0){
524  corrlumiAvg_h->SetBinError(bx,totalLumiByBXAvg_[bx]*correctionScaleFactors_[bx]/TMath::Sqrt(events_.at(bx)));
525  }else{
526  corrlumiAvg_h->SetBinError(bx,0.0);
527  }
528 
529  scaleFactorAvg_h->SetBinContent(bx,correctionScaleFactors_[bx]);
530  lumiAvg_h->SetBinContent(bx,totalLumiByBXAvg_[bx]);
531  }
532 
533  //Writing the corrections to SQL lite file for db
540 
541 
542  if( poolDbService.isAvailable() ){
543  poolDbService->writeOne<LumiCorrections>(pccCorrections,thisIOV,"LumiCorrectionsRcd");
544  }
545  else{
546  throw std::runtime_error("PoolDBService required.");
547  }
548 
549  delete pccCorrections;
550 
551  histoFile->cd();
552  corrlumiAvg_h->Write();
553  scaleFactorAvg_h->Write();
554  lumiAvg_h->Write();
555 
556  delete corrlumiAvg_h;
557  delete scaleFactorAvg_h;
558  delete lumiAvg_h;
559 
560  type1FracHist->Fill(type1Frac);
563 
564 
565  type1FracGraph->SetPoint(iBlock,thisIOV+approxLumiBlockSize_/2.0,type1Frac);
571 
572  edm::LogInfo("INFO") <<"iBlock type1Frac mean_type1_residual mean_type2_residual mean_type1_residual_unc mean_type2_residual_unc "<<iBlock<<" "<<type1Frac<<" "<<mean_type1_residual<<" "<<mean_type2_residual<<" "<<mean_type1_residual_unc<<" "<<mean_type2_residual_unc;
573 
574  type1Frac=0.0;
578  mean_type2_residual_unc=0;
579 
580  iBlock++;
581 
582  resetBlock();
583  }
584  histoFile->cd();
585  type1FracHist->Write();
586  type1resHist->Write();
587  type2resHist->Write();
588 
589  delete type1FracHist;
590  delete type1resHist;
591  delete type2resHist;
592 
593  delete[] histname1;
594  delete[] histname2;
595  delete[] histname3;
596 
598  delete lumiInfoMapIterator->second;
599  }
600  for(unsigned int lumiSection=0; lumiSection<lumiSections.size(); lumiSection++){
601  thisLS=lumiSections[lumiSection];
602  delete lumiInfoMapPerLS[thisLS];
603  }
604  lumiInfoMap.clear();
605  lumiInfoMapPerLS.clear();
606  lumiSections.clear();
607  lumiInfoCounter.clear();
608 }
std::vector< float > events_
void calculateCorrections(std::vector< float >, std::vector< float > &, float &)
RunID const & id() const
Definition: RunBase.h:39
RunNumber_t run() const
Definition: RunBase.h:40
void setType1Fraction(float type1frac)
RunNumber_t run() const
Definition: RunID.h:39
std::map< std::pair< unsigned int, unsigned int >, LumiInfo * >::iterator lumiInfoMapIterator
void setOverallCorrection(float overallCorrection)
static const unsigned int numBX
Definition: LumiConstants.h:9
unsigned int iBlock
unsigned int LuminosityBlockNumber_t
unsigned long long Time_t
Definition: Time.h:16
std::vector< float > totalLumiByBXAvg_
TGraphErrors * type1resGraph
std::map< std::pair< unsigned int, unsigned int >, LumiInfo * > lumiInfoMap
LumiCorrections * pccCorrections
bool isAvailable() const
Definition: Service.h:46
void writeOne(T *payload, Time_t time, const std::string &recordName, bool withlogging=false)
std::vector< float > correctionScaleFactors_
void setType2Residual(float type2res)
std::vector< unsigned int > lumiSections
TGraphErrors * type2resGraph
std::vector< float > rawlumiBX_
std::vector< float > totalLumiByBX_
edm::Service< cond::service::PoolDBOutputService > poolDbService
TGraphErrors * type1FracGraph
unsigned int approxLumiBlockSize_
void setType1Residual(float type1res)
std::map< std::pair< unsigned int, unsigned int >, unsigned int > lumiInfoCounter
std::vector< float > errOnLumiByBX_
void setCorrectionsBX(std::vector< float > &correctBX)
unsigned int thisLS
std::map< unsigned int, LumiInfo * > lumiInfoMapPerLS
void CorrPCCProducer::estimateType1Frac ( std::vector< float >  uncorrPCCPerBX,
float &  type1Frac 
)
private

Definition at line 194 of file CorrPCCProducer.cc.

References correctionTemplate_, evaluateCorrectionResiduals(), mps_fire::i, gen::k, mean_type1_residual, and LumiConstants::numBX.

Referenced by calculateCorrections().

194  {
195  std::vector<float> corrected_tmp_;
196  for(size_t i=0; i<uncorrPCCPerBX.size(); i++){
197  corrected_tmp_.push_back(uncorrPCCPerBX.at(i));
198  }
199 
200 
201  //Apply initial type 1 correction
202  for(size_t k=0;k<LumiConstants::numBX-1; k++){
203  float bin_k = corrected_tmp_.at(k);
204  corrected_tmp_.at(k+1)=corrected_tmp_.at(k+1)-type1Frac*bin_k;
205  }
206 
207 
208  //Apply type 2 correction
209  for(size_t i=0; i<LumiConstants::numBX-1; i++){
210  for(size_t j=i+1; j<i+LumiConstants::numBX-1; j++){
211  float bin_i = corrected_tmp_.at(i);
212  if (j<LumiConstants::numBX){
213  corrected_tmp_.at(j)=corrected_tmp_.at(j)-bin_i*correctionTemplate_.at(j-i);
214  } else {
215  corrected_tmp_.at(j-LumiConstants::numBX) = corrected_tmp_.at(j-LumiConstants::numBX)-bin_i*correctionTemplate_.at(j-i);
216  }
217  }
218  }
219 
220  //Apply additional iteration for type 1 correction
221  evaluateCorrectionResiduals(corrected_tmp_);
223 }
static const unsigned int numBX
Definition: LumiConstants.h:9
void evaluateCorrectionResiduals(std::vector< float >)
std::vector< float > correctionTemplate_
int k[5][pyjets_maxn]
void CorrPCCProducer::evaluateCorrectionResiduals ( std::vector< float >  corrected_tmp_)
private

Definition at line 227 of file CorrPCCProducer.cc.

References getMaximum(), histoFile, diffTreeTool::index, mean_type1_residual, mean_type1_residual_unc, mean_type2_residual, mean_type2_residual_unc, nTrain, LumiConstants::numBX, and electronIdCutBased_cfi::threshold.

Referenced by calculateCorrections(), and estimateType1Frac().

227  {
228  float lumiMax = getMaximum(corrected_tmp_);
229  float threshold = lumiMax*0.2;
230 
231  mean_type1_residual = 0;
232  mean_type2_residual = 0;
233  nTrain = 0;
234  float lumi=0;
235  std::vector<float> afterGlow;
236  TH1F type1("type1","",1000,-0.5,0.5);
237  TH1F type2("type2","",1000,-0.5,0.5);
238  for(size_t ibx=2; ibx<LumiConstants::numBX-5; ibx++){
239  lumi = corrected_tmp_.at(ibx);
240  afterGlow.clear();
241  afterGlow.push_back(corrected_tmp_.at(ibx+1));
242  afterGlow.push_back(corrected_tmp_.at(ibx+2));
243 
244  //Where type 1 and type 2 residuals are computed
245  if(lumi>threshold && afterGlow[0]<threshold && afterGlow[1]<threshold){
246  for(int index=3; index<6; index++){
247  float thisAfterGlow=corrected_tmp_.at(ibx+index);
248  if(thisAfterGlow<threshold){
249  afterGlow.push_back(thisAfterGlow);
250  } else {
251  break;
252  }
253  }
254  float thisType1=0;
255  float thisType2=0;
256  if(afterGlow.size()>1){
257  int nAfter=0;
258  for(unsigned int index=1; index<afterGlow.size(); index++){
259  thisType2+=afterGlow[index];
260  type2.Fill(afterGlow[index]/lumi);
261  nAfter++;
262  }
263  thisType2/=nAfter;
264  }
265  thisType1=(afterGlow[0]-thisType2)/lumi;
266 
267  type1.Fill(thisType1);
268  nTrain+=1;
269  }
270  }
271 
272  mean_type1_residual=type1.GetMean();//Calculate the mean value of the type 1 residual
273  mean_type2_residual=type2.GetMean();//Calculate the mean value of the type 2 residual
274  mean_type1_residual_unc=type1.GetMeanError();
275  mean_type2_residual_unc=type2.GetMeanError();
276 
277  histoFile->cd();
278  type1.Write();
279  type2.Write();
280 }
unsigned int nTrain
static const unsigned int numBX
Definition: LumiConstants.h:9
float getMaximum(std::vector< float >)
float CorrPCCProducer::getMaximum ( std::vector< float >  lumi_vector)
private

Definition at line 181 of file CorrPCCProducer.cc.

References mps_fire::i.

Referenced by calculateCorrections(), and evaluateCorrectionResiduals().

181  {
182  float max_lumi=0;
183  for(size_t i=0; i<lumi_vector.size(); i++){
184  if(lumi_vector.at(i)>max_lumi) max_lumi = lumi_vector.at(i);
185  }
186  return max_lumi;
187 }
void CorrPCCProducer::makeCorrectionTemplate ( )
private

Definition at line 173 of file CorrPCCProducer.cc.

References correctionTemplate_, JetChargeProducer_cfi::exp, LumiConstants::numBX, type2_a_, and type2_b_.

Referenced by CorrPCCProducer().

173  {
174  for(unsigned int bx=1;bx<LumiConstants::numBX;bx++){
175  correctionTemplate_.at(bx)=type2_a_*exp(-(float(bx)-1)*type2_b_);
176  }
177 }
static const unsigned int numBX
Definition: LumiConstants.h:9
std::vector< float > correctionTemplate_
void CorrPCCProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
finalprivate
void CorrPCCProducer::resetBlock ( )
private

Definition at line 611 of file CorrPCCProducer.cc.

References correctionScaleFactors_, events_, LumiConstants::numBX, totalLumiByBX_, and totalLumiByBXAvg_.

Referenced by CorrPCCProducer(), and endRunProduce().

611  {
612  for(unsigned int bx=0;bx<LumiConstants::numBX;bx++){
613  totalLumiByBX_[bx]=0;
614  totalLumiByBXAvg_[bx]=0;
615  events_[bx]=0;
616  correctionScaleFactors_[bx]=1.0;
617  }
618 }
std::vector< float > events_
static const unsigned int numBX
Definition: LumiConstants.h:9
std::vector< float > totalLumiByBXAvg_
std::vector< float > correctionScaleFactors_
std::vector< float > totalLumiByBX_

Member Data Documentation

unsigned int CorrPCCProducer::approxLumiBlockSize_
private

Definition at line 111 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), and endRunProduce().

std::vector<float> CorrPCCProducer::correctionScaleFactors_
private

Definition at line 80 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), endRunProduce(), and resetBlock().

std::vector<float> CorrPCCProducer::correctionTemplate_
private
TH1F* CorrPCCProducer::corrlumiAvg_h
private

Definition at line 92 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

unsigned int CorrPCCProducer::countLumi_
private

Definition at line 110 of file CorrPCCProducer.cc.

Referenced by beginLuminosityBlock(), and CorrPCCProducer().

std::vector<float> CorrPCCProducer::errOnLumiByBX_
private

Definition at line 75 of file CorrPCCProducer.cc.

Referenced by endLuminosityBlock(), and endRunProduce().

std::vector<float> CorrPCCProducer::events_
private

Definition at line 78 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), endLuminosityBlock(), endRunProduce(), and resetBlock().

TFile* CorrPCCProducer::histoFile
private
TList* CorrPCCProducer::hlist
private

Definition at line 101 of file CorrPCCProducer.cc.

unsigned int CorrPCCProducer::iBlock =0
private

Definition at line 83 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

TH1F* CorrPCCProducer::lumiAvg_h
private

Definition at line 94 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

std::map<std::pair<unsigned int,unsigned int>, unsigned int> CorrPCCProducer::lumiInfoCounter
private

Definition at line 90 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

std::map<std::pair<unsigned int,unsigned int>, LumiInfo*> CorrPCCProducer::lumiInfoMap
private

Definition at line 89 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

std::map<std::pair<unsigned int,unsigned int>, LumiInfo*>::iterator CorrPCCProducer::lumiInfoMapIterator
private

Definition at line 88 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

std::map<unsigned int, LumiInfo*> CorrPCCProducer::lumiInfoMapPerLS
private

Definition at line 86 of file CorrPCCProducer.cc.

Referenced by endLuminosityBlock(), and endRunProduce().

edm::EDGetTokenT<LumiInfo> CorrPCCProducer::lumiInfoToken
private

Definition at line 70 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), and endLuminosityBlock().

std::vector<unsigned int> CorrPCCProducer::lumiSections
private

Definition at line 87 of file CorrPCCProducer.cc.

Referenced by endLuminosityBlock(), and endRunProduce().

float CorrPCCProducer::mean_type1_residual
private
float CorrPCCProducer::mean_type1_residual_unc
private

Definition at line 107 of file CorrPCCProducer.cc.

Referenced by endRunProduce(), and evaluateCorrectionResiduals().

float CorrPCCProducer::mean_type2_residual
private

Definition at line 106 of file CorrPCCProducer.cc.

Referenced by endRunProduce(), and evaluateCorrectionResiduals().

float CorrPCCProducer::mean_type2_residual_unc
private

Definition at line 108 of file CorrPCCProducer.cc.

Referenced by endRunProduce(), and evaluateCorrectionResiduals().

unsigned int CorrPCCProducer::minimumNumberOfEvents
private

Definition at line 84 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), and endLuminosityBlock().

unsigned int CorrPCCProducer::nTrain
private

Definition at line 109 of file CorrPCCProducer.cc.

Referenced by evaluateCorrectionResiduals().

float CorrPCCProducer::overallCorrection_
private

Definition at line 81 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

LumiCorrections* CorrPCCProducer::pccCorrections
private

Definition at line 118 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

std::string CorrPCCProducer::pccSrc_
private

Definition at line 71 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer().

edm::Service<cond::service::PoolDBOutputService> CorrPCCProducer::poolDbService
private

Definition at line 120 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

std::string CorrPCCProducer::prodInst_
private

Definition at line 72 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer().

std::vector<float> CorrPCCProducer::rawlumiBX_
private

Definition at line 74 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

TH1F* CorrPCCProducer::scaleFactorAvg_h
private

Definition at line 93 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

unsigned int CorrPCCProducer::thisLS
private

Definition at line 112 of file CorrPCCProducer.cc.

Referenced by endLuminosityBlock(), and endRunProduce().

std::vector<float> CorrPCCProducer::totalLumiByBX_
private

Definition at line 76 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), endLuminosityBlock(), endRunProduce(), and resetBlock().

std::vector<float> CorrPCCProducer::totalLumiByBXAvg_
private

Definition at line 77 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), endRunProduce(), and resetBlock().

float CorrPCCProducer::type1Frac
private

Definition at line 104 of file CorrPCCProducer.cc.

Referenced by calculateCorrections(), and endRunProduce().

TGraphErrors* CorrPCCProducer::type1FracGraph
private

Definition at line 98 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), endJob(), and endRunProduce().

TH1F* CorrPCCProducer::type1FracHist
private

Definition at line 95 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

TGraphErrors* CorrPCCProducer::type1resGraph
private

Definition at line 99 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), endJob(), and endRunProduce().

TH1F* CorrPCCProducer::type1resHist
private

Definition at line 96 of file CorrPCCProducer.cc.

Referenced by endRunProduce().

double CorrPCCProducer::type2_a_
private

Definition at line 114 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), and makeCorrectionTemplate().

double CorrPCCProducer::type2_b_
private

Definition at line 115 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), and makeCorrectionTemplate().

TGraphErrors* CorrPCCProducer::type2resGraph
private

Definition at line 100 of file CorrPCCProducer.cc.

Referenced by CorrPCCProducer(), endJob(), and endRunProduce().

TH1F* CorrPCCProducer::type2resHist
private

Definition at line 97 of file CorrPCCProducer.cc.

Referenced by endRunProduce().