CMS 3D CMS Logo

Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes

GetLumi Class Reference

#include <GetLumi.h>

List of all members.

Public Types

enum  SQRT_S { SQRT_S_7TeV, SQRT_S_8TeV }

Public Member Functions

double convert2PU (double, double)
double convert2PU (double, int)
 GetLumi (const edm::InputTag &, double)
 GetLumi (const edm::ParameterSet &)
double getRawValue (const edm::Event &)
double getRawValue (edm::LuminosityBlock const &, edm::EventSetup const &)
double getValue (edm::LuminosityBlock const &, edm::EventSetup const &)
double getValue (const edm::Event &)
virtual ~GetLumi ()

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)

Public Attributes

edm::InputTag lumiInputTag_
double lumiScale_

Static Public Attributes

static constexpr double FREQ_ORBIT = 11246.
static constexpr double INELASTIC_XSEC_7TeV = 68.0
static constexpr double INELASTIC_XSEC_8TeV = 69.3
static const unsigned int NUM_BX = 3564
static constexpr double SECONDS_PER_LS = double(0x40000)/double(FREQ_ORBIT)

Detailed Description

Definition at line 16 of file GetLumi.h.


Member Enumeration Documentation

Enumerator:
SQRT_S_7TeV 
SQRT_S_8TeV 

Definition at line 20 of file GetLumi.h.


Constructor & Destructor Documentation

GetLumi::GetLumi ( const edm::ParameterSet iConfig)

Definition at line 15 of file GetLumi.cc.

  : lumiInputTag_ ( iConfig.getParameter<edm::InputTag>("lumi")  )
  , lumiScale_    ( iConfig.getParameter<double>("lumiScale")    )
{
}
GetLumi::GetLumi ( const edm::InputTag lumiInputTag,
double  lumiScale 
)

Definition at line 21 of file GetLumi.cc.

  : lumiInputTag_ ( lumiInputTag )
  , lumiScale_    ( lumiScale    )
{
}
GetLumi::~GetLumi ( ) [virtual]

Definition at line 27 of file GetLumi.cc.

{
}

Member Function Documentation

double GetLumi::convert2PU ( double  instLumi,
double  inelastic_xSec = GetLumi::INELASTIC_XSEC_8TeV 
)

Definition at line 94 of file GetLumi.cc.

References FREQ_ORBIT.

Referenced by convert2PU().

{

  // from https://cmswbm.web.cern.ch/cmswbm/images/pileup.png
  return instLumi*inelastic_xSec/FREQ_ORBIT;
}
double GetLumi::convert2PU ( double  instLumi,
int  sqrt_s = GetLumi::SQRT_S_8TeV 
)

Definition at line 102 of file GetLumi.cc.

References convert2PU(), INELASTIC_XSEC_7TeV, INELASTIC_XSEC_8TeV, SQRT_S_7TeV, and SQRT_S_8TeV.

{
  
  double inelastic_xSec = 0.;
  
  switch(sqrt_s) {
  case GetLumi::SQRT_S_7TeV :
    inelastic_xSec = GetLumi::INELASTIC_XSEC_7TeV;
    break;
  case GetLumi::SQRT_S_8TeV :
    inelastic_xSec = GetLumi::INELASTIC_XSEC_8TeV;
    break;
  default :
    break;
  }
  
  return convert2PU(instLumi,inelastic_xSec);

}
static void GetLumi::fillDescriptions ( edm::ConfigurationDescriptions descriptions) [static]
double GetLumi::getRawValue ( const edm::Event iEvent)

Definition at line 32 of file GetLumi.cc.

References edm::EventBase::bunchCrossing(), edm::Event::getLuminosityBlock(), edm::HandleBase::isValid(), LumiDetails::kOCC1, fjr2json::lumi, and lumiInputTag_.

Referenced by getValue().

{

  // taken from 
  // DPGAnalysis/SiStripTools/src/DigiLumiCorrHistogramMaker.cc
  // the scale factor 6.37 should follow the lumi prescriptions
  edm::Handle<LumiDetails> lumi;
  iEvent.getLuminosityBlock().getByLabel(lumiInputTag_,lumi);

  double bxlumi = -1.;
  if(lumi->isValid()) {
    bxlumi = lumi->lumiValue(LumiDetails::kOCC1,iEvent.bunchCrossing());
  }

  return bxlumi;

}
double GetLumi::getRawValue ( edm::LuminosityBlock const &  lumiBlock,
edm::EventSetup const &  eSetup 
)

Definition at line 59 of file GetLumi.cc.

References gather_cfg::cout, edm::LuminosityBlock::getByLabel(), edm::HandleBase::isValid(), fjr2json::lumi, and lumiInputTag_.

{
  
  double lumi = -1.;
  double intDelLumi = -1.;

  //  size_t LS = lumiBlock.luminosityBlockAuxiliary().luminosityBlock();
  // accumulate HF data at every LS as it is closed. 
  // note: lumi unit from DIPLumiSummary and Detail is microbarns
  edm::Handle<LumiSummary> lumiSummary_;
  lumiBlock.getByLabel(lumiInputTag_, lumiSummary_);
  if(lumiSummary_->isValid()){
    lumi = lumiSummary_->avgInsDelLumi();
    intDelLumi = lumiSummary_->intgDelLumi();
    std::cout << "Luminosity in this Lumi Section " << lumi << " --> " << intDelLumi << std::endl;
  } else {
    std::cout << "No valid data found!" << std::endl;
  }


  return lumi;

}
double GetLumi::getValue ( const edm::Event iEvent)

Definition at line 52 of file GetLumi.cc.

References getRawValue(), and lumiScale_.

Referenced by VertexMonitor::analyze(), TrackingMonitor::analyze(), and LogMessageMonitor::analyze().

{
  //    bxlumi = lumi->lumiValue(LumiDetails::kOCC1,iEvent.bunchCrossing())*6.37;
  return getRawValue(iEvent)*lumiScale_;
}
double GetLumi::getValue ( edm::LuminosityBlock const &  lumiBlock,
edm::EventSetup const &  eSetup 
)

Definition at line 86 of file GetLumi.cc.

References getRawValue(), and lumiScale_.

{
  return getRawValue(lumiBlock,eSetup)*lumiScale_;
}

Member Data Documentation

constexpr double GetLumi::FREQ_ORBIT = 11246. [static]

Definition at line 26 of file GetLumi.h.

Referenced by convert2PU().

constexpr double GetLumi::INELASTIC_XSEC_7TeV = 68.0 [static]

Definition at line 29 of file GetLumi.h.

Referenced by convert2PU().

constexpr double GetLumi::INELASTIC_XSEC_8TeV = 69.3 [static]

Definition at line 30 of file GetLumi.h.

Referenced by convert2PU().

Definition at line 47 of file GetLumi.h.

Referenced by getRawValue().

Definition at line 48 of file GetLumi.h.

Referenced by getValue().

const unsigned int GetLumi::NUM_BX = 3564 [static]

Definition at line 25 of file GetLumi.h.

constexpr double GetLumi::SECONDS_PER_LS = double(0x40000)/double(FREQ_ORBIT) [static]

Definition at line 27 of file GetLumi.h.