CMS 3D CMS Logo

GetLumi.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author: Mia Tosi,40 3-B32,+41227671609
5  */
6 
7 #include <iostream>
9 
13 
15  : lumiInputTag_(iConfig.getParameter<edm::InputTag>("lumi")),
16  lumiScale_(iConfig.getParameter<double>("lumiScale")) {}
17 
19  : lumiInputTag_(lumiInputTag), lumiScale_(lumiScale) {}
20 
24 }
25 
27  : GetLumi(lumiInputTag, lumiScale) {
30 }
31 
33 
35  // taken from
36  // DPGAnalysis/SiStripTools/src/DigiLumiCorrHistogramMaker.cc
37  // the scale factor 6.37 should follow the lumi prescriptions
40 
41  double bxlumi = -1.;
42  if (lumi->isValid()) {
43  bxlumi = lumi->lumiValue(LumiDetails::kOCC1, iEvent.bunchCrossing());
44  }
45 
46  return bxlumi;
47 }
48 
50  // bxlumi = lumi->lumiValue(LumiDetails::kOCC1,iEvent.bunchCrossing())*6.37;
51  return getRawValue(iEvent) * lumiScale_;
52 }
53 
54 double GetLumi::getRawValue(edm::LuminosityBlock const& lumiBlock, edm::EventSetup const& eSetup) {
55  double lumi = -1.;
56  double intDelLumi = -1.;
57 
58  // size_t LS = lumiBlock.luminosityBlockAuxiliary().luminosityBlock();
59  // accumulate HF data at every LS as it is closed.
60  // note: lumi unit from DIPLumiSummary and Detail is microbarns
61  edm::Handle<LumiSummary> lumiSummary_;
62  lumiBlock.getByToken(lumiSummaryToken_, lumiSummary_);
63  if (lumiSummary_->isValid()) {
64  lumi = lumiSummary_->avgInsDelLumi();
65  intDelLumi = lumiSummary_->intgDelLumi();
66  std::cout << "Luminosity in this Lumi Section " << lumi << " --> " << intDelLumi << std::endl;
67  } else {
68  std::cout << "No valid data found!" << std::endl;
69  }
70 
71  return lumi;
72 }
73 
74 double GetLumi::getValue(edm::LuminosityBlock const& lumiBlock, edm::EventSetup const& eSetup) {
75  return getRawValue(lumiBlock, eSetup) * lumiScale_;
76 }
77 
79  double inelastic_xSec = GetLumi::INELASTIC_XSEC_8TeV) // inelastic_xSec in mb
80 {
81  // from https://cmswbm.web.cern.ch/cmswbm/images/pileup.png
82  return instLumi * inelastic_xSec / FREQ_ORBIT;
83 }
84 
85 double GetLumi::convert2PU(double instLumi, int sqrt_s = GetLumi::SQRT_S_8TeV) {
86  double inelastic_xSec = 0.;
87 
88  switch (sqrt_s) {
90  inelastic_xSec = GetLumi::INELASTIC_XSEC_7TeV;
91  break;
93  inelastic_xSec = GetLumi::INELASTIC_XSEC_8TeV;
94  break;
95  default:
96  break;
97  }
98 
99  return convert2PU(instLumi, inelastic_xSec);
100 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
static double INELASTIC_XSEC_8TeV
Definition: GetLumi.h:30
bool isValid() const
Definition: LumiDetails.cc:44
static double INELASTIC_XSEC_7TeV
Definition: GetLumi.h:29
int bunchCrossing() const
Definition: EventBase.h:64
edm::EDGetTokenT< LumiDetails > lumiDetailsToken_
Definition: GetLumi.h:65
float lumiValue(AlgoType algo, unsigned int bx) const
Definition: LumiDetails.cc:72
virtual ~GetLumi()
Definition: GetLumi.cc:32
double getValue(const edm::Event &)
Definition: GetLumi.cc:49
int iEvent
Definition: GenABIO.cc:224
float intgDelLumi() const
Definition: LumiSummary.cc:16
edm::EDGetTokenT< LumiSummary > lumiSummaryToken_
Definition: GetLumi.h:66
LuminosityBlock const & getLuminosityBlock() const
Definition: Event.h:98
double getRawValue(const edm::Event &)
Definition: GetLumi.cc:34
float avgInsDelLumi() const
Definition: LumiSummary.cc:8
GetLumi(const edm::ParameterSet &)
Definition: GetLumi.cc:14
double convert2PU(double, double)
Definition: GetLumi.cc:78
HLT enums.
edm::InputTag lumiInputTag_
Definition: GetLumi.h:62
bool isValid() const
Definition: LumiSummary.cc:50
static double FREQ_ORBIT
Definition: GetLumi.h:26
double lumiScale_
Definition: GetLumi.h:63