CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes
HGCalConcentratorSuperTriggerCellImpl Class Reference

#include <HGCalConcentratorSuperTriggerCellImpl.h>

Classes

class  SuperTriggerCell
 

Public Member Functions

 HGCalConcentratorSuperTriggerCellImpl (const edm::ParameterSet &conf)
 
void select (const std::vector< l1t::HGCalTriggerCell > &trigCellVecInput, std::vector< l1t::HGCalTriggerCell > &trigCellVecOutput)
 
void setGeometry (const HGCalTriggerGeometryBase *const geom)
 

Private Types

enum  EnergyDivisionType { superTriggerCell, oneBitFraction, equalShare }
 

Private Member Functions

void assignSuperTriggerCellEnergyAndPosition (l1t::HGCalTriggerCell &c, const SuperTriggerCell &stc) const
 
void createAllTriggerCells (std::unordered_map< unsigned, SuperTriggerCell > &STCs, std::vector< l1t::HGCalTriggerCell > &trigCellVecOutput) const
 
uint32_t getCompressedSTCEnergy (const SuperTriggerCell &stc) const
 
float getTriggerCellOneBitFraction (float tcPt, float sumPt) const
 

Private Attributes

HGCalTriggerCellCalibration calibrationEE_
 
HGCalTriggerCellCalibration calibrationHEsc_
 
HGCalTriggerCellCalibration calibrationHEsi_
 
HGCalTriggerCellCalibration calibrationNose_
 
std::vector< unsigned > coarsenTriggerCells_
 
HGCalCoarseTriggerCellMapping coarseTCmapping_
 
EnergyDivisionType energyDivisionType_
 
bool fixedDataSizePerHGCROC_
 
double oneBitFractionHighValue_
 
double oneBitFractionLowValue_
 
double oneBitFractionThreshold_
 
HGCalCoarseTriggerCellMapping superTCmapping_
 
HGCalTriggerTools triggerTools_
 
HGCalVFECompressionImpl vfeCompression_
 

Static Private Attributes

static constexpr int kHighDensityThickness_ = 0
 
static constexpr int kOddNumberMask_ = 1
 
static constexpr int kTriggerCellsForDivision_ = 4
 

Detailed Description

Definition at line 13 of file HGCalConcentratorSuperTriggerCellImpl.h.

Member Enumeration Documentation

◆ EnergyDivisionType

Constructor & Destructor Documentation

◆ HGCalConcentratorSuperTriggerCellImpl()

HGCalConcentratorSuperTriggerCellImpl::HGCalConcentratorSuperTriggerCellImpl ( const edm::ParameterSet conf)

Definition at line 3 of file HGCalConcentratorSuperTriggerCellImpl.cc.

References energyDivisionType_, equalShare, edm::ParameterSet::getParameter(), oneBitFraction, oneBitFractionHighValue_, oneBitFractionLowValue_, oneBitFractionThreshold_, AlCaHLTBitMon_QueryRunRegistry::string, and superTriggerCell.

4  : fixedDataSizePerHGCROC_(conf.getParameter<bool>("fixedDataSizePerHGCROC")),
5  coarsenTriggerCells_(conf.getParameter<std::vector<unsigned>>("coarsenTriggerCells")),
6  coarseTCmapping_(conf.getParameter<std::vector<unsigned>>("ctcSize")),
7  superTCmapping_(conf.getParameter<std::vector<unsigned>>("stcSize")),
8  calibrationEE_(conf.getParameterSet("superTCCalibration_ee")),
9  calibrationHEsi_(conf.getParameterSet("superTCCalibration_hesi")),
10  calibrationHEsc_(conf.getParameterSet("superTCCalibration_hesc")),
11  calibrationNose_(conf.getParameterSet("superTCCalibration_nose")),
12  vfeCompression_(conf.getParameterSet("superTCCompression")) {
13  std::string energyType(conf.getParameter<string>("type_energy_division"));
14 
15  if (energyType == "superTriggerCell") {
17  } else if (energyType == "oneBitFraction") {
19 
20  oneBitFractionThreshold_ = conf.getParameter<double>("oneBitFractionThreshold");
21  oneBitFractionLowValue_ = conf.getParameter<double>("oneBitFractionLowValue");
22  oneBitFractionHighValue_ = conf.getParameter<double>("oneBitFractionHighValue");
23 
24  } else if (energyType == "equalShare") {
26 
27  } else {
29  }
30 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
ParameterSet const & getParameterSet(std::string const &) const

Member Function Documentation

◆ assignSuperTriggerCellEnergyAndPosition()

void HGCalConcentratorSuperTriggerCellImpl::assignSuperTriggerCellEnergyAndPosition ( l1t::HGCalTriggerCell c,
const SuperTriggerCell stc 
) const
private

Definition at line 92 of file HGCalConcentratorSuperTriggerCellImpl.cc.

References HltBtagPostValidation_cff::c, HGCalTriggerCellCalibration::calibrateInGeV(), calibrationEE_, calibrationHEsc_, calibrationHEsi_, calibrationNose_, coarsenTriggerCells_, coarseTCmapping_, bTagMiniDQMDeepCSV::denominator, energyDivisionType_, equalShare, Exception, fixedDataSizePerHGCROC_, DivergingColor::frac, HGCalCoarseTriggerCellMapping::getCoarseTriggerCellId(), HGCalCoarseTriggerCellMapping::getCoarseTriggerCellPosition(), getCompressedSTCEnergy(), HGCalCoarseTriggerCellMapping::getConstituentTriggerCells(), HGCalConcentratorSuperTriggerCellImpl::SuperTriggerCell::getFractionSum(), HGCalConcentratorSuperTriggerCellImpl::SuperTriggerCell::getMaxId(), HGCalConcentratorSuperTriggerCellImpl::SuperTriggerCell::getSTCId(), HGCalTriggerTools::getSubDetectorType(), HGCalConcentratorSuperTriggerCellImpl::SuperTriggerCell::getSumPt(), HGCalTriggerTools::getTCPosition(), HGCalConcentratorSuperTriggerCellImpl::SuperTriggerCell::getTCpt(), getTriggerCellOneBitFraction(), HGCalTriggerTools::isEm(), HGCalTriggerTools::isNose(), HGCalTriggerTools::isSilicon(), kHighDensityThickness_, kTriggerCellsForDivision_, oneBitFraction, point, findQualityFiles::size, superTCmapping_, superTriggerCell, Calorimetry_cff::thickness, HGCalTriggerTools::thicknessIndex(), and triggerTools_.

Referenced by createAllTriggerCells().

93  {
94  //Compress and recalibrate STC energy
95  uint32_t compressed_value = getCompressedSTCEnergy(stc);
96 
98  int thickness = triggerTools_.thicknessIndex(c.detId());
99 
100  bool isSilicon = triggerTools_.isSilicon(c.detId());
101  bool isEM = triggerTools_.isEm(c.detId());
102  bool isNose = triggerTools_.isNose(c.detId());
103 
107  } else {
108  point = triggerTools_.getTCPosition(c.detId());
109  }
110  c.setPosition(point);
111 
112  math::PtEtaPhiMLorentzVector p4(c.pt(), point.eta(), point.phi(), 0.);
113  c.setP4(p4);
114 
116  if (c.detId() == stc.getMaxId()) {
117  c.setHwPt(compressed_value);
118  } else {
119  throw cms::Exception("NonMaxIdSuperTriggerCell")
120  << "Trigger Cell with detId not equal to the maximum of the superTriggerCell found";
121  }
122  } else if (energyDivisionType_ == equalShare) {
123  double coarseTriggerCellSize =
124  coarsenTriggerCells_[subdet]
125  ? double(
127  .size())
128  : 1.;
129 
130  double denominator =
132  ? double(kTriggerCellsForDivision_)
133  : double(superTCmapping_.getConstituentTriggerCells(stc.getSTCId()).size()) / coarseTriggerCellSize;
134 
135  c.setHwPt(std::round(compressed_value / denominator));
136 
137  } else if (energyDivisionType_ == oneBitFraction) {
138  double frac = 0;
139 
140  if (c.detId() != stc.getMaxId()) {
141  frac = getTriggerCellOneBitFraction(stc.getTCpt(c.detId()), stc.getSumPt());
142  } else {
143  frac = 1 - stc.getFractionSum();
144  }
145 
146  c.setHwPt(std::round(compressed_value * frac));
147  }
148  // calibration
149  if (isNose) {
151  } else if (isSilicon) {
152  if (isEM) {
154  } else {
156  }
157  } else {
159  }
160 }
size
Write out results.
std::vector< uint32_t > getConstituentTriggerCells(uint32_t ctcId) const
uint32_t getCompressedSTCEnergy(const SuperTriggerCell &stc) const
bool isSilicon(const DetId &) const
int thicknessIndex(const DetId &) const
void calibrateInGeV(l1t::HGCalTriggerCell &) const
uint32_t getCoarseTriggerCellId(uint32_t detid) const
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
GlobalPoint getTCPosition(const DetId &id) const
SubDetectorType getSubDetectorType(const DetId &id) const
GlobalPoint getCoarseTriggerCellPosition(uint32_t ctcId) const
bool isNose(const DetId &) const
bool isEm(const DetId &) const
float getTriggerCellOneBitFraction(float tcPt, float sumPt) const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5

◆ createAllTriggerCells()

void HGCalConcentratorSuperTriggerCellImpl::createAllTriggerCells ( std::unordered_map< unsigned, SuperTriggerCell > &  STCs,
std::vector< l1t::HGCalTriggerCell > &  trigCellVecOutput 
) const
private

Definition at line 43 of file HGCalConcentratorSuperTriggerCellImpl.cc.

References assignSuperTriggerCellEnergyAndPosition(), coarsenTriggerCells_, coarseTCmapping_, energyDivisionType_, Exception, fixedDataSizePerHGCROC_, HGCalCoarseTriggerCellMapping::getCoarseTriggerCellId(), HGCalCoarseTriggerCellMapping::getConstituentTriggerCells(), HGCalCoarseTriggerCellMapping::getRepresentativeDetId(), HGCalTriggerTools::getSubDetectorType(), getTriggerCellOneBitFraction(), HGCalTriggerTools::getTriggerGeometry(), kHighDensityThickness_, oneBitFraction, alignCSCRings::s, l1t::HGCalTriggerCell::setDetId(), superTCmapping_, superTriggerCell, Calorimetry_cff::thickness, HGCalTriggerTools::thicknessIndex(), triggerTools_, and HGCalTriggerGeometryBase::validTriggerCell().

Referenced by select().

44  {
45  for (auto& s : STCs) {
46  std::vector<uint32_t> output_ids = superTCmapping_.getConstituentTriggerCells(s.second.getSTCId());
47  if (output_ids.empty())
48  continue;
49 
51  int thickness = (!output_ids.empty() ? triggerTools_.thicknessIndex(output_ids.at(0)) : 0);
52 
53  for (const auto& id : output_ids) {
56  continue;
57  }
58 
60  continue;
61  }
62 
63  l1t::HGCalTriggerCell triggerCell;
64  triggerCell.setDetId(id);
65  if (energyDivisionType_ == superTriggerCell && id != s.second.getMaxId()) {
66  continue;
67  }
68 
69  DetId tc_Id(id);
70 
71  if (superTCmapping_.getCoarseTriggerCellId(id) != s.second.getSTCId()) {
72  throw cms::Exception("NonExistingCoarseTC")
73  << "The coarse trigger cell correponsing to the nominal trigger cell does not exist";
74  }
75  trigCellVecOutput.push_back(triggerCell);
76  if (energyDivisionType_ == oneBitFraction) { //Get the 1 bit fractions
77 
78  if (id != s.second.getMaxId()) {
79  float tc_fraction = getTriggerCellOneBitFraction(s.second.getTCpt(id), s.second.getSumPt());
80  s.second.addToFractionSum(tc_fraction);
81  }
82  }
83  }
84  }
85  // assign energy
86  for (l1t::HGCalTriggerCell& tc : trigCellVecOutput) {
87  const auto& stc = STCs[superTCmapping_.getCoarseTriggerCellId(tc.detId())];
89  }
90 }
std::vector< uint32_t > getConstituentTriggerCells(uint32_t ctcId) const
uint32_t getRepresentativeDetId(uint32_t tcid) const
int thicknessIndex(const DetId &) const
const HGCalTriggerGeometryBase * getTriggerGeometry() const
uint32_t getCoarseTriggerCellId(uint32_t detid) const
virtual bool validTriggerCell(const unsigned trigger_cell_id) const =0
void setDetId(uint32_t detid)
SubDetectorType getSubDetectorType(const DetId &id) const
void assignSuperTriggerCellEnergyAndPosition(l1t::HGCalTriggerCell &c, const SuperTriggerCell &stc) const
Definition: DetId.h:17
float getTriggerCellOneBitFraction(float tcPt, float sumPt) const

◆ getCompressedSTCEnergy()

uint32_t HGCalConcentratorSuperTriggerCellImpl::getCompressedSTCEnergy ( const SuperTriggerCell stc) const
private

Definition at line 32 of file HGCalConcentratorSuperTriggerCellImpl.cc.

References HGCalVFECompressionImpl::compressSingle(), HGCalConcentratorSuperTriggerCellImpl::SuperTriggerCell::getSumHwPt(), SiStripPI::max, and vfeCompression_.

Referenced by assignSuperTriggerCellEnergyAndPosition().

32  {
33  uint32_t code(0);
34  uint64_t compressed_value(0);
35  vfeCompression_.compressSingle(stc.getSumHwPt(), code, compressed_value);
36 
37  if (compressed_value > std::numeric_limits<uint32_t>::max())
38  edm::LogWarning("CompressedValueDowncasting") << "Compressed value cannot fit into 32-bit word. Downcasting.";
39 
40  return static_cast<uint32_t>(compressed_value);
41 }
void compressSingle(const uint64_t value, uint32_t &compressedCode, uint64_t &compressedValue) const
unsigned long long uint64_t
Definition: Time.h:13

◆ getTriggerCellOneBitFraction()

float HGCalConcentratorSuperTriggerCellImpl::getTriggerCellOneBitFraction ( float  tcPt,
float  sumPt 
) const
private

◆ select()

void HGCalConcentratorSuperTriggerCellImpl::select ( const std::vector< l1t::HGCalTriggerCell > &  trigCellVecInput,
std::vector< l1t::HGCalTriggerCell > &  trigCellVecOutput 
)

Definition at line 174 of file HGCalConcentratorSuperTriggerCellImpl.cc.

References createAllTriggerCells(), HGCalCoarseTriggerCellMapping::getCoarseTriggerCellId(), and superTCmapping_.

175  {
176  std::unordered_map<unsigned, SuperTriggerCell> STCs;
177  // first pass, fill the "coarse" trigger cells
178  for (const l1t::HGCalTriggerCell& tc : trigCellVecInput) {
179  uint32_t stcid = superTCmapping_.getCoarseTriggerCellId(tc.detId());
180  STCs[stcid].add(tc, stcid);
181  }
182 
183  createAllTriggerCells(STCs, trigCellVecOutput);
184 }
void createAllTriggerCells(std::unordered_map< unsigned, SuperTriggerCell > &STCs, std::vector< l1t::HGCalTriggerCell > &trigCellVecOutput) const
uint32_t getCoarseTriggerCellId(uint32_t detid) const

◆ setGeometry()

void HGCalConcentratorSuperTriggerCellImpl::setGeometry ( const HGCalTriggerGeometryBase *const  geom)
inline

Definition at line 19 of file HGCalConcentratorSuperTriggerCellImpl.h.

References calibrationEE_, calibrationHEsc_, calibrationHEsi_, calibrationNose_, coarseTCmapping_, relativeConstraints::geom, HGCalTriggerCellCalibration::setGeometry(), HGCalCoarseTriggerCellMapping::setGeometry(), HGCalTriggerTools::setGeometry(), superTCmapping_, and triggerTools_.

19  {
27  }
void setGeometry(const HGCalTriggerGeometryBase *const geom)
void setGeometry(const HGCalTriggerGeometryBase *const)
void setGeometry(const HGCalTriggerGeometryBase *const geom)

Member Data Documentation

◆ calibrationEE_

HGCalTriggerCellCalibration HGCalConcentratorSuperTriggerCellImpl::calibrationEE_
private

◆ calibrationHEsc_

HGCalTriggerCellCalibration HGCalConcentratorSuperTriggerCellImpl::calibrationHEsc_
private

◆ calibrationHEsi_

HGCalTriggerCellCalibration HGCalConcentratorSuperTriggerCellImpl::calibrationHEsi_
private

◆ calibrationNose_

HGCalTriggerCellCalibration HGCalConcentratorSuperTriggerCellImpl::calibrationNose_
private

◆ coarsenTriggerCells_

std::vector<unsigned> HGCalConcentratorSuperTriggerCellImpl::coarsenTriggerCells_
private

◆ coarseTCmapping_

HGCalCoarseTriggerCellMapping HGCalConcentratorSuperTriggerCellImpl::coarseTCmapping_
private

◆ energyDivisionType_

EnergyDivisionType HGCalConcentratorSuperTriggerCellImpl::energyDivisionType_
private

◆ fixedDataSizePerHGCROC_

bool HGCalConcentratorSuperTriggerCellImpl::fixedDataSizePerHGCROC_
private

◆ kHighDensityThickness_

constexpr int HGCalConcentratorSuperTriggerCellImpl::kHighDensityThickness_ = 0
staticprivate

◆ kOddNumberMask_

constexpr int HGCalConcentratorSuperTriggerCellImpl::kOddNumberMask_ = 1
staticprivate

Definition at line 38 of file HGCalConcentratorSuperTriggerCellImpl.h.

◆ kTriggerCellsForDivision_

constexpr int HGCalConcentratorSuperTriggerCellImpl::kTriggerCellsForDivision_ = 4
staticprivate

◆ oneBitFractionHighValue_

double HGCalConcentratorSuperTriggerCellImpl::oneBitFractionHighValue_
private

◆ oneBitFractionLowValue_

double HGCalConcentratorSuperTriggerCellImpl::oneBitFractionLowValue_
private

◆ oneBitFractionThreshold_

double HGCalConcentratorSuperTriggerCellImpl::oneBitFractionThreshold_
private

◆ superTCmapping_

HGCalCoarseTriggerCellMapping HGCalConcentratorSuperTriggerCellImpl::superTCmapping_
private

◆ triggerTools_

HGCalTriggerTools HGCalConcentratorSuperTriggerCellImpl::triggerTools_
private

◆ vfeCompression_

HGCalVFECompressionImpl HGCalConcentratorSuperTriggerCellImpl::vfeCompression_
private

Definition at line 58 of file HGCalConcentratorSuperTriggerCellImpl.h.

Referenced by getCompressedSTCEnergy().