CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
HGCalVFESummationImpl Class Reference

#include <HGCalVFESummationImpl.h>

Public Member Functions

void eventSetup (const edm::EventSetup &es)
 
 HGCalVFESummationImpl (const edm::ParameterSet &conf)
 
void triggerCellSums (const HGCalTriggerGeometryBase &, const std::vector< std::pair< DetId, uint32_t > > &, std::unordered_map< uint32_t, uint32_t > &payload)
 

Private Attributes

double lsb_scintillator_MIP_
 
double lsb_silicon_fC_
 
std::vector< double > thickness_corrections_
 
double threshold_scintillator_
 
std::vector< double > thresholds_silicon_
 
HGCalTriggerTools triggerTools_
 

Detailed Description

Definition at line 14 of file HGCalVFESummationImpl.h.

Constructor & Destructor Documentation

HGCalVFESummationImpl::HGCalVFESummationImpl ( const edm::ParameterSet conf)

Definition at line 3 of file HGCalVFESummationImpl.cc.

References Exception, edm::ParameterSet::getParameter(), hgcalDigitizer_cfi::noise, thickness_corrections_, MessageLogger_cff::threshold, threshold_scintillator_, thresholds_silicon_, and HcalDetIdTransform::transform().

4  : thickness_corrections_(conf.getParameter<std::vector<double>>("ThicknessCorrections")),
5  lsb_silicon_fC_(conf.getParameter<double>("siliconCellLSB_fC")),
6  lsb_scintillator_MIP_(conf.getParameter<double>("scintillatorCellLSB_MIP")) {
7  const unsigned nThickness = 3;
8  if (thickness_corrections_.size() != nThickness) {
9  throw cms::Exception("Configuration")
10  << thickness_corrections_.size() << " thickness corrections are given instead of " << nThickness
11  << " (the number of sensor thicknesses)";
12  }
14  conf.getParameter<edm::ParameterSet>("noiseSilicon").getParameter<std::vector<double>>("values");
15  if (thresholds_silicon_.size() != nThickness) {
16  throw cms::Exception("Configuration") << thresholds_silicon_.size() << " silicon thresholds are given instead of "
17  << nThickness << " (the number of sensor thicknesses)";
18  }
19  threshold_scintillator_ = conf.getParameter<edm::ParameterSet>("noiseScintillator").getParameter<double>("noise_MIP");
20  const auto threshold = conf.getParameter<double>("noiseThreshold");
23  return noise * threshold;
24  });
26 }
T getParameter(std::string const &) const
std::vector< double > thresholds_silicon_
std::vector< double > thickness_corrections_
unsigned transform(const HcalDetId &id, unsigned transformCode)

Member Function Documentation

void HGCalVFESummationImpl::eventSetup ( const edm::EventSetup es)
inline

Definition at line 18 of file HGCalVFESummationImpl.h.

References HGCalTriggerTools::eventSetup(), jets_cff::payload, triggerCellSums(), and triggerTools_.

void eventSetup(const edm::EventSetup &)
HGCalTriggerTools triggerTools_
void HGCalVFESummationImpl::triggerCellSums ( const HGCalTriggerGeometryBase ,
const std::vector< std::pair< DetId, uint32_t > > &  ,
std::unordered_map< uint32_t, uint32_t > &  payload 
)

Definition at line 28 of file HGCalVFESummationImpl.cc.

References amptDefault_cfi::frame, HGCalTriggerGeometryBase::getTriggerCellFromCell(), HGCalTriggerTools::isScintillator(), HGCalTriggerTools::isSilicon(), lsb_scintillator_MIP_, lsb_silicon_fC_, Calorimetry_cff::thickness, thickness_corrections_, HGCalTriggerTools::thicknessIndex(), MessageLogger_cff::threshold, threshold_scintillator_, thresholds_silicon_, triggerTools_, and relativeConstraints::value.

Referenced by eventSetup().

30  {
31  if (linearized_dataframes.empty())
32  return;
33  // sum energies in trigger cells
34  for (const auto& frame : linearized_dataframes) {
35  DetId cellid(frame.first);
36  uint32_t value = frame.second;
37 
38  // Apply noise threshold before summing into trigger cells
39  if (triggerTools_.isSilicon(cellid)) {
41  double threshold = thresholds_silicon_.at(thickness);
42  value = (value * lsb_silicon_fC_ > threshold ? value : 0);
43  } else if (triggerTools_.isScintillator(cellid)) {
44  value = (value * lsb_scintillator_MIP_ > threshold_scintillator_ ? value : 0);
45  }
46  if (value == 0)
47  continue;
48 
49  // find trigger cell associated to cell
50  uint32_t tcid = geometry.getTriggerCellFromCell(cellid);
51  payload.emplace(tcid, 0); // do nothing if key exists already
52 
53  // equalize value among cell thicknesses for Silicon parts
54  if (triggerTools_.isSilicon(cellid)) {
55  int thickness = triggerTools_.thicknessIndex(cellid);
56  double thickness_correction = thickness_corrections_.at(thickness);
57  value = (double)value * thickness_correction;
58  }
59 
60  // sums energy for the same trigger cell id
61  payload[tcid] += value; // 32 bits integer should be largely enough
62  }
63 }
bool isScintillator(const DetId &id) const
HGCalTriggerTools triggerTools_
std::vector< double > thresholds_silicon_
std::vector< double > thickness_corrections_
Definition: value.py:1
Definition: DetId.h:17
int thicknessIndex(const DetId &, bool tc=false) const
bool isSilicon(const DetId &) const

Member Data Documentation

double HGCalVFESummationImpl::lsb_scintillator_MIP_
private

Definition at line 26 of file HGCalVFESummationImpl.h.

Referenced by triggerCellSums().

double HGCalVFESummationImpl::lsb_silicon_fC_
private

Definition at line 25 of file HGCalVFESummationImpl.h.

Referenced by triggerCellSums().

std::vector<double> HGCalVFESummationImpl::thickness_corrections_
private

Definition at line 24 of file HGCalVFESummationImpl.h.

Referenced by HGCalVFESummationImpl(), and triggerCellSums().

double HGCalVFESummationImpl::threshold_scintillator_
private

Definition at line 28 of file HGCalVFESummationImpl.h.

Referenced by HGCalVFESummationImpl(), and triggerCellSums().

std::vector<double> HGCalVFESummationImpl::thresholds_silicon_
private

Definition at line 27 of file HGCalVFESummationImpl.h.

Referenced by HGCalVFESummationImpl(), and triggerCellSums().

HGCalTriggerTools HGCalVFESummationImpl::triggerTools_
private

Definition at line 29 of file HGCalVFESummationImpl.h.

Referenced by eventSetup(), and triggerCellSums().