CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCalConcentratorSuperTriggerCellImpl.cc
Go to the documentation of this file.
2 
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 }
31 
33  uint32_t code(0);
34  uint32_t compressed_value(0);
35  vfeCompression_.compressSingle(stc.getSumHwPt(), code, compressed_value);
36  return compressed_value;
37 }
38 
40  std::unordered_map<unsigned, SuperTriggerCell>& STCs, std::vector<l1t::HGCalTriggerCell>& trigCellVecOutput) const {
41  for (auto& s : STCs) {
42  std::vector<uint32_t> output_ids = superTCmapping_.getConstituentTriggerCells(s.second.getSTCId());
43  if (output_ids.empty())
44  continue;
45 
47  int thickness = (!output_ids.empty() ? triggerTools_.thicknessIndex(output_ids.at(0)) : 0);
48 
49  for (const auto& id : output_ids) {
50  if (((fixedDataSizePerHGCROC_ && thickness > kHighDensityThickness_) || coarsenTriggerCells_[subdet]) &&
52  continue;
53  }
54 
56  continue;
57  }
58 
59  l1t::HGCalTriggerCell triggerCell;
60  triggerCell.setDetId(id);
61  if (energyDivisionType_ == superTriggerCell && id != s.second.getMaxId()) {
62  continue;
63  }
64 
65  DetId tc_Id(id);
66 
67  if (superTCmapping_.getCoarseTriggerCellId(id) != s.second.getSTCId()) {
68  throw cms::Exception("NonExistingCoarseTC")
69  << "The coarse trigger cell correponsing to the nominal trigger cell does not exist";
70  }
71  trigCellVecOutput.push_back(triggerCell);
72  if (energyDivisionType_ == oneBitFraction) { //Get the 1 bit fractions
73 
74  if (id != s.second.getMaxId()) {
75  float tc_fraction = getTriggerCellOneBitFraction(s.second.getTCpt(id), s.second.getSumPt());
76  s.second.addToFractionSum(tc_fraction);
77  }
78  }
79  }
80  }
81  // assign energy
82  for (l1t::HGCalTriggerCell& tc : trigCellVecOutput) {
83  const auto& stc = STCs[superTCmapping_.getCoarseTriggerCellId(tc.detId())];
85  }
86 }
87 
89  const SuperTriggerCell& stc) const {
90  //Compress and recalibrate STC energy
91  uint32_t compressed_value = getCompressedSTCEnergy(stc);
92 
95 
96  bool isSilicon = triggerTools_.isSilicon(c.detId());
97  bool isEM = triggerTools_.isEm(c.detId());
98  bool isNose = triggerTools_.isNose(c.detId());
99 
101  if ((fixedDataSizePerHGCROC_ && thickness > kHighDensityThickness_) || coarsenTriggerCells_[subdet]) {
103  } else {
104  point = triggerTools_.getTCPosition(c.detId());
105  }
106  c.setPosition(point);
107 
108  math::PtEtaPhiMLorentzVector p4(c.pt(), point.eta(), point.phi(), 0.);
109  c.setP4(p4);
110 
112  if (c.detId() == stc.getMaxId()) {
113  c.setHwPt(compressed_value);
114  } else {
115  throw cms::Exception("NonMaxIdSuperTriggerCell")
116  << "Trigger Cell with detId not equal to the maximum of the superTriggerCell found";
117  }
118  } else if (energyDivisionType_ == equalShare) {
119  double coarseTriggerCellSize =
120  coarsenTriggerCells_[subdet]
121  ? double(
123  .size())
124  : 1.;
125 
126  double denominator =
128  ? double(kTriggerCellsForDivision_)
129  : double(superTCmapping_.getConstituentTriggerCells(stc.getSTCId()).size()) / coarseTriggerCellSize;
130 
131  c.setHwPt(std::round(compressed_value / denominator));
132 
133  } else if (energyDivisionType_ == oneBitFraction) {
134  double frac = 0;
135 
136  if (c.detId() != stc.getMaxId()) {
137  frac = getTriggerCellOneBitFraction(stc.getTCpt(c.detId()), stc.getSumPt());
138  } else {
139  frac = 1 - stc.getFractionSum();
140  }
141 
142  c.setHwPt(std::round(compressed_value * frac));
143  }
144  // calibration
145  if (isNose) {
147  } else if (isSilicon) {
148  if (isEM) {
150  } else {
152  }
153  } else {
155  }
156 }
157 
159  double f = tcPt / sumPt;
160  double frac = 0;
161  if (f < oneBitFractionThreshold_) {
163  } else {
165  }
166 
167  return frac;
168 }
169 
170 void HGCalConcentratorSuperTriggerCellImpl::select(const std::vector<l1t::HGCalTriggerCell>& trigCellVecInput,
171  std::vector<l1t::HGCalTriggerCell>& trigCellVecOutput) {
172  std::unordered_map<unsigned, SuperTriggerCell> STCs;
173  // first pass, fill the "coarse" trigger cells
174  for (const l1t::HGCalTriggerCell& tc : trigCellVecInput) {
175  uint32_t stcid = superTCmapping_.getCoarseTriggerCellId(tc.detId());
176  STCs[stcid].add(tc, stcid);
177  }
178 
179  createAllTriggerCells(STCs, trigCellVecOutput);
180 }
int thicknessIndex(const DetId &) const
const edm::EventSetup & c
float getTriggerCellOneBitFraction(float tcPt, float sumPt) const
double pt() const final
transverse momentum
ParameterSet const & getParameterSet(ParameterSetID const &id)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
const HGCalTriggerGeometryBase * getTriggerGeometry() const
uint32_t getCoarseTriggerCellId(uint32_t detid) const
void select(const std::vector< l1t::HGCalTriggerCell > &trigCellVecInput, std::vector< l1t::HGCalTriggerCell > &trigCellVecOutput)
HGCalConcentratorSuperTriggerCellImpl(const edm::ParameterSet &conf)
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
uint32_t detId() const
virtual bool validTriggerCell(const unsigned trigger_cell_id) const =0
list denominator
Definition: cuy.py:485
void setDetId(uint32_t detid)
bool isNose(const DetId &) const
Definition: DetId.h:17
bool isSilicon(const DetId &) const
void calibrateInGeV(l1t::HGCalTriggerCell &) const
uint32_t getRepresentativeDetId(uint32_t tcid) const
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
SubDetectorType getSubDetectorType(const DetId &id) const
T eta() const
Definition: PV3DBase.h:73
bool isEm(const DetId &) const
void setHwPt(int pt)
Definition: L1Candidate.h:28
std::vector< uint32_t > getConstituentTriggerCells(uint32_t ctcId) const
uint32_t getCompressedSTCEnergy(const SuperTriggerCell &stc) const
GlobalPoint getTCPosition(const DetId &id) const
void compressSingle(const uint32_t value, uint32_t &compressedCode, uint32_t &compressedValue) const
GlobalPoint getCoarseTriggerCellPosition(uint32_t ctcId) const
void createAllTriggerCells(std::unordered_map< unsigned, SuperTriggerCell > &STCs, std::vector< l1t::HGCalTriggerCell > &trigCellVecOutput) const
void setPosition(const GlobalPoint &position)
void setP4(const LorentzVector &p4) final
set 4-momentum
tuple size
Write out results.
*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
void assignSuperTriggerCellEnergyAndPosition(l1t::HGCalTriggerCell &c, const SuperTriggerCell &stc) const