CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
HGCalConcentratorCoarsenerImpl Class Reference

#include <HGCalConcentratorCoarsenerImpl.h>

Classes

struct  CoarseTC
 

Public Member Functions

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

Private Member Functions

void assignCoarseTriggerCellEnergy (l1t::HGCalTriggerCell &c, const CoarseTC &ctc) const
 
void updateCoarseTriggerCellMaps (const l1t::HGCalTriggerCell &tc, uint32_t ctcid)
 

Private Attributes

HGCalTriggerCellCalibration calibration_
 
HGCalCoarseTriggerCellMapping coarseTCmapping_
 
std::unordered_map< uint32_t,
CoarseTC
coarseTCs_
 
bool fixedDataSizePerHGCROC_
 
HGCalTriggerTools triggerTools_
 
HGCalVFECompressionImpl vfeCompression_
 

Static Private Attributes

static constexpr int kHighDensityThickness_ = 0
 

Detailed Description

Definition at line 10 of file HGCalConcentratorCoarsenerImpl.h.

Constructor & Destructor Documentation

HGCalConcentratorCoarsenerImpl::HGCalConcentratorCoarsenerImpl ( const edm::ParameterSet conf)

Definition at line 3 of file HGCalConcentratorCoarsenerImpl.cc.

4  : fixedDataSizePerHGCROC_(conf.getParameter<bool>("fixedDataSizePerHGCROC")),
5  coarseTCmapping_(conf.getParameter<std::vector<unsigned>>("ctcSize")),
6  calibration_(conf.getParameterSet("superTCCalibration")),
7  vfeCompression_(conf.getParameterSet("coarseTCCompression")) {}
ParameterSet const & getParameterSet(std::string const &) const
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HGCalCoarseTriggerCellMapping coarseTCmapping_

Member Function Documentation

void HGCalConcentratorCoarsenerImpl::assignCoarseTriggerCellEnergy ( l1t::HGCalTriggerCell c,
const CoarseTC ctc 
) const
private

Definition at line 22 of file HGCalConcentratorCoarsenerImpl.cc.

References HGCalTriggerCellCalibration::calibrateInGeV(), calibration_, HGCalVFECompressionImpl::compressSingle(), l1t::L1Candidate::setHwPt(), HGCalConcentratorCoarsenerImpl::CoarseTC::sumHwPt, and vfeCompression_.

Referenced by coarsen().

23  {
24  //Compress and recalibrate CTC energy
25  uint32_t code(0);
26  uint32_t compressed_value(0);
27  vfeCompression_.compressSingle(ctc.sumHwPt, code, compressed_value);
28 
29  tc.setHwPt(compressed_value);
31 }
void calibrateInGeV(l1t::HGCalTriggerCell &) const
void compressSingle(const uint32_t value, uint32_t &compressedCode, uint32_t &compressedValue) const
void HGCalConcentratorCoarsenerImpl::coarsen ( const std::vector< l1t::HGCalTriggerCell > &  trigCellVecInput,
std::vector< l1t::HGCalTriggerCell > &  trigCellVecOutput 
)

Definition at line 33 of file HGCalConcentratorCoarsenerImpl.cc.

References assignCoarseTriggerCellEnergy(), coarseTCmapping_, coarseTCs_, PV3DBase< T, PVType, FrameType >::eta(), fixedDataSizePerHGCROC_, HGCalCoarseTriggerCellMapping::getCoarseTriggerCellId(), HGCalCoarseTriggerCellMapping::getCoarseTriggerCellPosition(), HGCalCoarseTriggerCellMapping::getRepresentativeDetId(), kHighDensityThickness_, PV3DBase< T, PVType, FrameType >::phi(), point, reco::LeafCandidate::pt(), l1t::HGCalTriggerCell::setDetId(), reco::LeafCandidate::setP4(), l1t::HGCalTriggerCell::setPosition(), TrackerMaterial_cfi::thickness, HGCalTriggerTools::thicknessIndex(), triggerTools_, and updateCoarseTriggerCellMaps().

34  {
35  coarseTCs_.clear();
36 
37  // first pass, fill the coarse trigger cell information
38  for (const l1t::HGCalTriggerCell& tc : trigCellVecInput) {
39  int thickness = triggerTools_.thicknessIndex(tc.detId());
40 
41  if (fixedDataSizePerHGCROC_ && thickness == kHighDensityThickness_) {
42  trigCellVecOutput.push_back(tc);
43  continue;
44  }
45 
46  uint32_t ctcid = coarseTCmapping_.getCoarseTriggerCellId(tc.detId());
47  updateCoarseTriggerCellMaps(tc, ctcid);
48  }
49 
50  for (const auto& ctc : coarseTCs_) {
51  l1t::HGCalTriggerCell triggerCell;
52 
53  uint32_t representativeId = coarseTCmapping_.getRepresentativeDetId(ctc.second.maxId);
54  triggerCell.setDetId(representativeId);
55 
57  math::PtEtaPhiMLorentzVector initial_p4(ctc.second.sumPt, point.eta(), point.phi(), 0);
58 
59  triggerCell.setP4(initial_p4);
60  assignCoarseTriggerCellEnergy(triggerCell, ctc.second);
61 
62  math::PtEtaPhiMLorentzVector p4(triggerCell.pt(), point.eta(), point.phi(), 0);
63  triggerCell.setPosition(point);
64  triggerCell.setP4(p4);
65  trigCellVecOutput.push_back(triggerCell);
66  }
67 }
int thicknessIndex(const DetId &) const
double pt() const final
transverse momentum
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
void updateCoarseTriggerCellMaps(const l1t::HGCalTriggerCell &tc, uint32_t ctcid)
uint32_t getCoarseTriggerCellId(uint32_t detid) const
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
void setDetId(uint32_t detid)
uint32_t getRepresentativeDetId(uint32_t tcid) const
T eta() const
Definition: PV3DBase.h:73
void assignCoarseTriggerCellEnergy(l1t::HGCalTriggerCell &c, const CoarseTC &ctc) const
std::unordered_map< uint32_t, CoarseTC > coarseTCs_
HGCalCoarseTriggerCellMapping coarseTCmapping_
GlobalPoint getCoarseTriggerCellPosition(uint32_t ctcId) const
void setPosition(const GlobalPoint &position)
void setP4(const LorentzVector &p4) final
set 4-momentum
*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 HGCalConcentratorCoarsenerImpl::setGeometry ( const HGCalTriggerGeometryBase *const  geom)
inline

Definition at line 16 of file HGCalConcentratorCoarsenerImpl.h.

References calibration_, coarseTCmapping_, HGCalTriggerCellCalibration::setGeometry(), HGCalCoarseTriggerCellMapping::setGeometry(), HGCalTriggerTools::setGeometry(), and triggerTools_.

16  {
20  }
void setGeometry(const HGCalTriggerGeometryBase *const geom)
void setGeometry(const HGCalTriggerGeometryBase *const)
HGCalCoarseTriggerCellMapping coarseTCmapping_
void setGeometry(const HGCalTriggerGeometryBase *const geom)
void HGCalConcentratorCoarsenerImpl::updateCoarseTriggerCellMaps ( const l1t::HGCalTriggerCell tc,
uint32_t  ctcid 
)
private

Definition at line 9 of file HGCalConcentratorCoarsenerImpl.cc.

References coarseTCs_, l1t::HGCalTriggerCell::detId(), l1t::L1Candidate::hwPt(), l1t::HGCalTriggerCell::mipPt(), and reco::LeafCandidate::pt().

Referenced by coarsen().

9  {
10  auto& ctc = coarseTCs_[ctcid];
11 
12  ctc.sumPt += tc.pt();
13  ctc.sumHwPt += tc.hwPt();
14  ctc.sumMipPt += tc.mipPt();
15 
16  if (tc.mipPt() > ctc.maxMipPt) {
17  ctc.maxId = tc.detId();
18  ctc.maxMipPt = tc.mipPt();
19  }
20 }
double pt() const final
transverse momentum
double mipPt() const
uint32_t detId() const
int hwPt() const
Definition: L1Candidate.h:35
std::unordered_map< uint32_t, CoarseTC > coarseTCs_

Member Data Documentation

HGCalTriggerCellCalibration HGCalConcentratorCoarsenerImpl::calibration_
private

Definition at line 28 of file HGCalConcentratorCoarsenerImpl.h.

Referenced by assignCoarseTriggerCellEnergy(), and setGeometry().

HGCalCoarseTriggerCellMapping HGCalConcentratorCoarsenerImpl::coarseTCmapping_
private

Definition at line 25 of file HGCalConcentratorCoarsenerImpl.h.

Referenced by coarsen(), and setGeometry().

std::unordered_map<uint32_t, CoarseTC> HGCalConcentratorCoarsenerImpl::coarseTCs_
private

Definition at line 39 of file HGCalConcentratorCoarsenerImpl.h.

Referenced by coarsen(), and updateCoarseTriggerCellMaps().

bool HGCalConcentratorCoarsenerImpl::fixedDataSizePerHGCROC_
private

Definition at line 24 of file HGCalConcentratorCoarsenerImpl.h.

Referenced by coarsen().

constexpr int HGCalConcentratorCoarsenerImpl::kHighDensityThickness_ = 0
staticprivate

Definition at line 26 of file HGCalConcentratorCoarsenerImpl.h.

Referenced by coarsen().

HGCalTriggerTools HGCalConcentratorCoarsenerImpl::triggerTools_
private

Definition at line 23 of file HGCalConcentratorCoarsenerImpl.h.

Referenced by coarsen(), and setGeometry().

HGCalVFECompressionImpl HGCalConcentratorCoarsenerImpl::vfeCompression_
private

Definition at line 29 of file HGCalConcentratorCoarsenerImpl.h.

Referenced by assignCoarseTriggerCellEnergy().