CMS 3D CMS Logo

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

#include <HGCalVFEProcessorSums.h>

Inheritance diagram for HGCalVFEProcessorSums:
HGCalProcessorBaseT< InputCollection, OutputCollection >

Public Member Functions

 HGCalVFEProcessorSums (const edm::ParameterSet &conf)
 
void run (const HGCalDigiCollection &digiColl, l1t::HGCalTriggerCellBxCollection &triggerCellColl) override
 
- Public Member Functions inherited from HGCalProcessorBaseT< InputCollection, OutputCollection >
 HGCalProcessorBaseT (const edm::ParameterSet &conf)
 
const std::string & name () const
 
virtual void run (const InputCollection &inputColl, OutputCollection &outColl)=0
 
virtual void setGeometry (const HGCalTriggerGeometryBase *const geom)
 
virtual ~HGCalProcessorBaseT ()
 

Private Attributes

std::unique_ptr< HGCalTriggerCellCalibrationcalibrationEE_
 
std::unique_ptr< HGCalTriggerCellCalibrationcalibrationHEsc_
 
std::unique_ptr< HGCalTriggerCellCalibrationcalibrationHEsi_
 
std::unique_ptr< HGCalTriggerCellCalibrationcalibrationNose_
 
HGCalTriggerTools triggerTools_
 
std::unique_ptr< HGCalVFECompressionImplvfeCompressionHDMImpl_
 
std::unique_ptr< HGCalVFECompressionImplvfeCompressionLDMImpl_
 
std::unique_ptr< HGCalVFELinearizationImplvfeLinearizationScImpl_
 
std::unique_ptr< HGCalVFELinearizationImplvfeLinearizationSiImpl_
 
std::unique_ptr< HGCalVFESummationImplvfeSummationImpl_
 

Additional Inherited Members

- Protected Member Functions inherited from HGCalProcessorBaseT< InputCollection, OutputCollection >
const HGCalTriggerGeometryBasegeometry () const
 

Detailed Description

Definition at line 12 of file HGCalVFEProcessorSums.h.

Constructor & Destructor Documentation

◆ HGCalVFEProcessorSums()

HGCalVFEProcessorSums::HGCalVFEProcessorSums ( const edm::ParameterSet conf)

Definition at line 5 of file HGCalVFEProcessorSums.cc.

References calibrationEE_, calibrationHEsc_, calibrationHEsi_, calibrationNose_, edm::ParameterSet::getParameter(), vfeCompressionHDMImpl_, vfeCompressionLDMImpl_, vfeLinearizationScImpl_, vfeLinearizationSiImpl_, and vfeSummationImpl_.

5  : HGCalVFEProcessorBase(conf) {
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 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< HGCalVFECompressionImpl > vfeCompressionHDMImpl_
HGCalProcessorBaseT< HGCalDigiCollection, l1t::HGCalTriggerCellBxCollection > HGCalVFEProcessorBase
std::unique_ptr< HGCalTriggerCellCalibration > calibrationHEsi_
std::unique_ptr< HGCalTriggerCellCalibration > calibrationHEsc_
std::unique_ptr< HGCalVFECompressionImpl > vfeCompressionLDMImpl_
std::unique_ptr< HGCalVFELinearizationImpl > vfeLinearizationScImpl_
std::unique_ptr< HGCalVFELinearizationImpl > vfeLinearizationSiImpl_
std::unique_ptr< HGCalTriggerCellCalibration > calibrationNose_
std::unique_ptr< HGCalVFESummationImpl > vfeSummationImpl_
std::unique_ptr< HGCalTriggerCellCalibration > calibrationEE_

Member Function Documentation

◆ run()

void HGCalVFEProcessorSums::run ( const HGCalDigiCollection digiColl,
l1t::HGCalTriggerCellBxCollection triggerCellColl 
)
override

Definition at line 28 of file HGCalVFEProcessorSums.cc.

References calibrationEE_, calibrationHEsc_, calibrationHEsi_, calibrationNose_, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), DigiToRawDM_cff::digiColl, HGCalTriggerGeometryBase::disconnectedModule(), HGCalProcessorBaseT< InputCollection, OutputCollection >::geometry(), HGCalTriggerGeometryBase::getModuleFromCell(), HGCalTriggerGeometryBase::getTriggerCellPosition(), HFNose, l1t::L1Candidate::hwPt(), mps_fire::i, HGCalTriggerTools::isEm(), HGCalTriggerTools::isNose(), HGCalTriggerTools::isSilicon(), WZElectronSkims53X_cff::max, point, BXVector< T >::push_back(), l1t::HGCalTriggerCell::setCompressedCharge(), HGCalTriggerTools::setGeometry(), reco::LeafCandidate::setP4(), l1t::HGCalTriggerCell::setPosition(), l1t::HGCalTriggerCell::setUncompressedCharge(), ppsPixelTopologyESSourceRun2_cfi::thickness, HGCalTriggerTools::thicknessIndex(), triggerTools_, vfeCompressionHDMImpl_, vfeCompressionLDMImpl_, vfeLinearizationScImpl_, vfeLinearizationSiImpl_, and vfeSummationImpl_.

29  {
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<uint64_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 
87  if (tc_compressed_value > std::numeric_limits<int>::max())
88  edm::LogWarning("CompressedValueDowncasting") << "Compressed value cannot fit into 32-bit word. Downcasting.";
89 
90  l1t::HGCalTriggerCell triggerCell(
91  reco::LeafCandidate::LorentzVector(), static_cast<int>(tc_compressed_value), 0, 0, 0, tc_id);
92 
93  if (tc_compressed_code > std::numeric_limits<uint32_t>::max())
94  edm::LogWarning("CompressedValueDowncasting") << "Compressed code cannot fit into 32-bit word. Downcasting.";
95 
96  triggerCell.setCompressedCharge(static_cast<uint32_t>(tc_compressed_code));
97  triggerCell.setUncompressedCharge(tc_value);
99 
100  // 'value' is hardware, so p4 is meaningless, except for eta and phi
101  math::PtEtaPhiMLorentzVector p4((double)tc_compressed_value / cosh(point.eta()), point.eta(), point.phi(), 0.);
102  triggerCell.setP4(p4);
103  triggerCell.setPosition(point);
104 
105  // calibration
106  if (triggerCell.hwPt() > 0) {
107  l1t::HGCalTriggerCell calibratedtriggercell(triggerCell);
108  if (isNose) {
109  calibrationNose_->calibrateInGeV(calibratedtriggercell);
110  } else if (isSilicon) {
111  if (isEM) {
112  calibrationEE_->calibrateInGeV(calibratedtriggercell);
113  } else {
114  calibrationHEsi_->calibrateInGeV(calibratedtriggercell);
115  }
116  } else {
117  calibrationHEsc_->calibrateInGeV(calibratedtriggercell);
118  }
119  triggerCellColl.push_back(0, calibratedtriggercell);
120  }
121  }
122  }
123 }
virtual GlobalPoint getTriggerCellPosition(const unsigned trigger_cell_det_id) const =0
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
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
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
std::unique_ptr< HGCalVFESummationImpl > vfeSummationImpl_
std::unique_ptr< HGCalTriggerCellCalibration > calibrationEE_
Log< level::Warning, false > LogWarning
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)

Member Data Documentation

◆ calibrationEE_

std::unique_ptr<HGCalTriggerCellCalibration> HGCalVFEProcessorSums::calibrationEE_
private

Definition at line 24 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().

◆ calibrationHEsc_

std::unique_ptr<HGCalTriggerCellCalibration> HGCalVFEProcessorSums::calibrationHEsc_
private

Definition at line 26 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().

◆ calibrationHEsi_

std::unique_ptr<HGCalTriggerCellCalibration> HGCalVFEProcessorSums::calibrationHEsi_
private

Definition at line 25 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().

◆ calibrationNose_

std::unique_ptr<HGCalTriggerCellCalibration> HGCalVFEProcessorSums::calibrationNose_
private

Definition at line 27 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().

◆ triggerTools_

HGCalTriggerTools HGCalVFEProcessorSums::triggerTools_
private

Definition at line 29 of file HGCalVFEProcessorSums.h.

Referenced by run().

◆ vfeCompressionHDMImpl_

std::unique_ptr<HGCalVFECompressionImpl> HGCalVFEProcessorSums::vfeCompressionHDMImpl_
private

Definition at line 23 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().

◆ vfeCompressionLDMImpl_

std::unique_ptr<HGCalVFECompressionImpl> HGCalVFEProcessorSums::vfeCompressionLDMImpl_
private

Definition at line 22 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().

◆ vfeLinearizationScImpl_

std::unique_ptr<HGCalVFELinearizationImpl> HGCalVFEProcessorSums::vfeLinearizationScImpl_
private

Definition at line 20 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().

◆ vfeLinearizationSiImpl_

std::unique_ptr<HGCalVFELinearizationImpl> HGCalVFEProcessorSums::vfeLinearizationSiImpl_
private

Definition at line 19 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().

◆ vfeSummationImpl_

std::unique_ptr<HGCalVFESummationImpl> HGCalVFEProcessorSums::vfeSummationImpl_
private

Definition at line 21 of file HGCalVFEProcessorSums.h.

Referenced by HGCalVFEProcessorSums(), and run().