CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 {
18 }
19 
20 GetLumi::GetLumi(const edm::InputTag& lumiInputTag, double lumiScale)
21  : lumiInputTag_ ( lumiInputTag )
22  , lumiScale_ ( lumiScale )
23 {
24 }
25 
27 {
30 }
31 
32 GetLumi::GetLumi(const edm::InputTag& lumiInputTag, double lumiScale, edm::ConsumesCollector& iC) : GetLumi(lumiInputTag,lumiScale)
33 {
36 }
37 
39 {
40 }
41 
42 double
44 {
45 
46  // taken from
47  // DPGAnalysis/SiStripTools/src/DigiLumiCorrHistogramMaker.cc
48  // the scale factor 6.37 should follow the lumi prescriptions
50  // iEvent.getLuminosityBlock().getByLabel(lumiInputTag_,lumi);
52 
53 
54  double bxlumi = -1.;
55  if(lumi->isValid()) {
56  bxlumi = lumi->lumiValue(LumiDetails::kOCC1,iEvent.bunchCrossing());
57  }
58 
59  return bxlumi;
60 
61 }
62 
63 
64 double
66 {
67  // bxlumi = lumi->lumiValue(LumiDetails::kOCC1,iEvent.bunchCrossing())*6.37;
68  return getRawValue(iEvent)*lumiScale_;
69 }
70 
71 double
73  edm::EventSetup const&eSetup)
74 {
75 
76  double lumi = -1.;
77  double intDelLumi = -1.;
78 
79  // size_t LS = lumiBlock.luminosityBlockAuxiliary().luminosityBlock();
80  // accumulate HF data at every LS as it is closed.
81  // note: lumi unit from DIPLumiSummary and Detail is microbarns
82  edm::Handle<LumiSummary> lumiSummary_;
83  // lumiBlock.getByLabel(lumiInputTag_, lumiSummary_);
84  lumiBlock.getByToken(lumiSummaryToken_, lumiSummary_);
85  if(lumiSummary_->isValid()){
86  lumi = lumiSummary_->avgInsDelLumi();
87  intDelLumi = lumiSummary_->intgDelLumi();
88  std::cout << "Luminosity in this Lumi Section " << lumi << " --> " << intDelLumi << std::endl;
89  } else {
90  std::cout << "No valid data found!" << std::endl;
91  }
92 
93 
94  return lumi;
95 
96 }
97 
98 
99 double
101  edm::EventSetup const&eSetup)
102 {
103  return getRawValue(lumiBlock,eSetup)*lumiScale_;
104 }
105 
106 
107 double
108 GetLumi::convert2PU(double instLumi, double inelastic_xSec = GetLumi::INELASTIC_XSEC_8TeV) // inelastic_xSec in mb
109 {
110 
111  // from https://cmswbm.web.cern.ch/cmswbm/images/pileup.png
112  return instLumi*inelastic_xSec/FREQ_ORBIT;
113 }
114 
115 double
116 GetLumi::convert2PU(double instLumi, int sqrt_s = GetLumi::SQRT_S_8TeV)
117 {
118 
119  double inelastic_xSec = 0.;
120 
121  switch(sqrt_s) {
122  case GetLumi::SQRT_S_7TeV :
123  inelastic_xSec = GetLumi::INELASTIC_XSEC_7TeV;
124  break;
125  case GetLumi::SQRT_S_8TeV :
126  inelastic_xSec = GetLumi::INELASTIC_XSEC_8TeV;
127  break;
128  default :
129  break;
130  }
131 
132  return convert2PU(instLumi,inelastic_xSec);
133 
134 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
tuple lumi
Definition: fjr2json.py:35
int bunchCrossing() const
Definition: EventBase.h:62
static constexpr double FREQ_ORBIT
Definition: GetLumi.h:31
edm::EDGetTokenT< LumiDetails > lumiDetailsToken_
Definition: GetLumi.h:57
virtual ~GetLumi()
Definition: GetLumi.cc:38
double getValue(const edm::Event &)
Definition: GetLumi.cc:65
int iEvent
Definition: GenABIO.cc:243
edm::EDGetTokenT< LumiSummary > lumiSummaryToken_
Definition: GetLumi.h:58
LuminosityBlock const & getLuminosityBlock() const
Definition: Event.h:80
bool isValid() const
Definition: HandleBase.h:76
static constexpr double INELASTIC_XSEC_7TeV
Definition: GetLumi.h:34
double getRawValue(const edm::Event &)
Definition: GetLumi.cc:43
static constexpr double INELASTIC_XSEC_8TeV
Definition: GetLumi.h:35
GetLumi(const edm::ParameterSet &)
Definition: GetLumi.cc:14
double convert2PU(double, double)
Definition: GetLumi.cc:108
tuple cout
Definition: gather_cfg.py:121
edm::InputTag lumiInputTag_
Definition: GetLumi.h:54
double lumiScale_
Definition: GetLumi.h:55