CMS 3D CMS Logo

HGCalVFEProcessorSums.cc
Go to the documentation of this file.
2 
4 
7  std::make_unique<HGCalVFELinearizationImpl>(conf.getParameter<edm::ParameterSet>("linearizationCfg_si"));
9  std::make_unique<HGCalVFELinearizationImpl>(conf.getParameter<edm::ParameterSet>("linearizationCfg_sc"));
10 
11  vfeSummationImpl_ = std::make_unique<HGCalVFESummationImpl>(conf.getParameter<edm::ParameterSet>("summationCfg"));
12 
14  std::make_unique<HGCalVFECompressionImpl>(conf.getParameter<edm::ParameterSet>("compressionCfg_ldm"));
16  std::make_unique<HGCalVFECompressionImpl>(conf.getParameter<edm::ParameterSet>("compressionCfg_hdm"));
17 
19  std::make_unique<HGCalTriggerCellCalibration>(conf.getParameter<edm::ParameterSet>("calibrationCfg_ee"));
21  std::make_unique<HGCalTriggerCellCalibration>(conf.getParameter<edm::ParameterSet>("calibrationCfg_hesi"));
23  std::make_unique<HGCalTriggerCellCalibration>(conf.getParameter<edm::ParameterSet>("calibrationCfg_hesc"));
25  std::make_unique<HGCalTriggerCellCalibration>(conf.getParameter<edm::ParameterSet>("calibrationCfg_nose"));
26 }
27 
29  l1t::HGCalTriggerCellBxCollection& triggerCellColl) {
30  vfeSummationImpl_->setGeometry(geometry());
31  calibrationEE_->setGeometry(geometry());
32  calibrationHEsi_->setGeometry(geometry());
33  calibrationHEsc_->setGeometry(geometry());
34  calibrationNose_->setGeometry(geometry());
36 
37  std::vector<HGCalDataFrame> dataframes;
38  std::vector<std::pair<DetId, uint32_t>> linearized_dataframes;
39  std::unordered_map<uint32_t, uint32_t> tc_payload;
40  std::unordered_map<uint32_t, std::array<uint32_t, 2>> tc_compressed_payload;
41 
42  // Remove disconnected modules and invalid cells
43  for (const auto& digiData : digiColl) {
44  if (!geometry()->validCell(digiData.id()))
45  continue;
46  uint32_t module = geometry()->getModuleFromCell(digiData.id());
47 
48  // no disconnected layer for HFNose
49  if (DetId(digiData.id()).subdetId() != ForwardSubdetector::HFNose) {
51  continue;
52  }
53 
54  dataframes.emplace_back(digiData.id());
55  for (int i = 0; i < digiData.size(); i++) {
56  dataframes.back().setSample(i, digiData.sample(i));
57  }
58  }
59  if (dataframes.empty())
60  return;
61 
62  constexpr int kHighDensityThickness = 0;
63  bool isSilicon = triggerTools_.isSilicon(dataframes[0].id());
64  bool isEM = triggerTools_.isEm(dataframes[0].id());
65  bool isNose = triggerTools_.isNose(dataframes[0].id());
66  int thickness = triggerTools_.thicknessIndex(dataframes[0].id());
67  // Linearization of ADC and TOT values to the same LSB
68  if (isSilicon) {
69  vfeLinearizationSiImpl_->linearize(dataframes, linearized_dataframes);
70  } else {
71  vfeLinearizationScImpl_->linearize(dataframes, linearized_dataframes);
72  }
73  // Sum of sensor cells into trigger cells
74  vfeSummationImpl_->triggerCellSums(linearized_dataframes, tc_payload);
75  // Compression of trigger cell charges to a floating point format
76  if (thickness == kHighDensityThickness) {
77  vfeCompressionHDMImpl_->compress(tc_payload, tc_compressed_payload);
78  } else {
79  vfeCompressionLDMImpl_->compress(tc_payload, tc_compressed_payload);
80  }
81 
82  // Transform map to trigger cell vector
83  for (const auto& [tc_id, tc_value] : tc_payload) {
84  if (tc_value > 0) {
85  const auto& [tc_compressed_code, tc_compressed_value] = tc_compressed_payload[tc_id];
86  l1t::HGCalTriggerCell triggerCell(reco::LeafCandidate::LorentzVector(), tc_compressed_value, 0, 0, 0, tc_id);
87  triggerCell.setCompressedCharge(tc_compressed_code);
88  triggerCell.setUncompressedCharge(tc_value);
90 
91  // 'value' is hardware, so p4 is meaningless, except for eta and phi
92  math::PtEtaPhiMLorentzVector p4((double)tc_compressed_value / cosh(point.eta()), point.eta(), point.phi(), 0.);
93  triggerCell.setP4(p4);
94  triggerCell.setPosition(point);
95 
96  // calibration
97  if (triggerCell.hwPt() > 0) {
98  l1t::HGCalTriggerCell calibratedtriggercell(triggerCell);
99  if (isNose) {
100  calibrationNose_->calibrateInGeV(calibratedtriggercell);
101  } else if (isSilicon) {
102  if (isEM) {
103  calibrationEE_->calibrateInGeV(calibratedtriggercell);
104  } else {
105  calibrationHEsi_->calibrateInGeV(calibratedtriggercell);
106  }
107  } else {
108  calibrationHEsc_->calibrateInGeV(calibratedtriggercell);
109  }
110  triggerCellColl.push_back(0, calibratedtriggercell);
111  }
112  }
113  }
114 }
virtual GlobalPoint getTriggerCellPosition(const unsigned trigger_cell_det_id) const =0
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
virtual unsigned getModuleFromCell(const unsigned cell_det_id) const =0
bool isSilicon(const DetId &) const
std::unique_ptr< HGCalVFECompressionImpl > vfeCompressionHDMImpl_
int thicknessIndex(const DetId &) const
HGCalVFEProcessorSums(const edm::ParameterSet &conf)
virtual bool disconnectedModule(const unsigned module_id) const =0
std::unique_ptr< HGCalTriggerCellCalibration > calibrationHEsi_
void setGeometry(const HGCalTriggerGeometryBase *const)
std::unique_ptr< HGCalTriggerCellCalibration > calibrationHEsc_
HGCalTriggerTools triggerTools_
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
std::unique_ptr< HGCalVFECompressionImpl > vfeCompressionLDMImpl_
const HGCalTriggerGeometryBase * geometry() const
void setCompressedCharge(uint32_t value)
bool isNose(const DetId &) const
bool isEm(const DetId &) const
std::unique_ptr< HGCalVFELinearizationImpl > vfeLinearizationScImpl_
std::unique_ptr< HGCalVFELinearizationImpl > vfeLinearizationSiImpl_
std::unique_ptr< HGCalTriggerCellCalibration > calibrationNose_
Definition: DetId.h:17
int hwPt() const
Definition: L1Candidate.h:35
std::unique_ptr< HGCalVFESummationImpl > vfeSummationImpl_
void setUncompressedCharge(uint32_t value)
std::unique_ptr< HGCalTriggerCellCalibration > calibrationEE_
#define DEFINE_EDM_PLUGIN(factory, type, name)
void setPosition(const GlobalPoint &position)
void setP4(const LorentzVector &p4) final
set 4-momentum
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
*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 push_back(int bx, T object)
void run(const HGCalDigiCollection &digiColl, l1t::HGCalTriggerCellBxCollection &triggerCellColl) override