CMS 3D CMS Logo

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

#include <HGCalTriggerCellCalibration.h>

Public Member Functions

void calibrateInGeV (l1t::HGCalTriggerCell &) const
 
void calibrateInMipT (l1t::HGCalTriggerCell &) const
 
void calibrateMipTinGeV (l1t::HGCalTriggerCell &) const
 
 HGCalTriggerCellCalibration (const edm::ParameterSet &conf)
 
void setGeometry (const HGCalTriggerGeometryBase *const geom)
 

Private Attributes

std::vector< double > chargeCollectionEfficiency_
 
std::vector< double > dEdX_weights_
 
std::vector< double > fCperMIP_
 
double lsb_
 
std::vector< double > thicknessCorrection_
 
HGCalTriggerTools triggerTools_
 

Detailed Description

Definition at line 9 of file HGCalTriggerCellCalibration.h.

Constructor & Destructor Documentation

◆ HGCalTriggerCellCalibration()

HGCalTriggerCellCalibration::HGCalTriggerCellCalibration ( const edm::ParameterSet conf)

Definition at line 5 of file HGCalTriggerCellCalibration.cc.

References chargeCollectionEfficiency_, l1tHGCalVFEProducer_cfi::fCperMIP, fCperMIP_, and thicknessCorrection_.

6  : lsb_(conf.getParameter<double>("lsb")),
7  fCperMIP_(conf.getParameter<std::vector<double>>("fCperMIP")),
8  chargeCollectionEfficiency_(conf.getParameter<edm::ParameterSet>("chargeCollectionEfficiency")
9  .getParameter<std::vector<double>>("values")),
10  thicknessCorrection_(conf.getParameter<std::vector<double>>("thicknessCorrection")),
11  dEdX_weights_(conf.getParameter<std::vector<double>>("dEdXweights")) {
12  for (const auto& fCperMIP : fCperMIP_) {
13  if (fCperMIP <= 0) {
14  edm::LogWarning("DivisionByZero") << "WARNING: zero or negative MIP->fC correction factor. It won't be "
15  "applied to correct trigger cell energies.";
16  }
17  }
18  for (const auto& cce : chargeCollectionEfficiency_) {
19  if (cce <= 0) {
20  edm::LogWarning("DivisionByZero") << "WARNING: zero or negative cell-thickness correction factor. It won't be "
21  "applied to correct trigger cell energies.";
22  }
23  }
24  for (const auto& thickCorr : thicknessCorrection_) {
25  if (thickCorr <= 0) {
26  edm::LogWarning("DivisionByZero") << "WARNING: zero or negative cell-thickness correction factor. It won't be "
27  "applied to correct trigger cell energies.";
28  }
29  }
30 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< double > chargeCollectionEfficiency_
std::vector< double > thicknessCorrection_
Log< level::Warning, false > LogWarning

Member Function Documentation

◆ calibrateInGeV()

void HGCalTriggerCellCalibration::calibrateInGeV ( l1t::HGCalTriggerCell trgCell) const

Definition at line 99 of file HGCalTriggerCellCalibration.cc.

References calibrateInMipT(), and calibrateMipTinGeV().

Referenced by HGCalConcentratorCoarsenerImpl::assignCoarseTriggerCellEnergy(), and HGCalConcentratorSuperTriggerCellImpl::assignSuperTriggerCellEnergyAndPosition().

99  {
100  /* calibrate from ADC count to transverse mip */
101  calibrateInMipT(trgCell);
102 
103  /* calibrate from mip count to GeV */
104  calibrateMipTinGeV(trgCell);
105 }
void calibrateInMipT(l1t::HGCalTriggerCell &) const
void calibrateMipTinGeV(l1t::HGCalTriggerCell &) const

◆ calibrateInMipT()

void HGCalTriggerCellCalibration::calibrateInMipT ( l1t::HGCalTriggerCell trgCell) const

Definition at line 32 of file HGCalTriggerCellCalibration.cc.

References CustomPhysics_cfi::amplitude, chargeCollectionEfficiency_, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), l1t::HGCalTriggerCell::detId(), reco::LeafCandidate::eta(), Exception, fCperMIP_, l1trig_cff::hwPt, l1t::L1Candidate::hwPt(), HGCalTriggerTools::isSilicon(), lsb_, l1t::HGCalTriggerCell::setMipPt(), Calorimetry_cff::thickness, HGCalTriggerTools::thicknessIndex(), and triggerTools_.

Referenced by calibrateInGeV().

32  {
33  DetId trgdetid(trgCell.detId());
34  bool isSilicon = triggerTools_.isSilicon(trgdetid);
35  constexpr int kScintillatorIndex = 0;
36  unsigned thickness = isSilicon ? triggerTools_.thicknessIndex(trgdetid) : kScintillatorIndex;
37  if (thickness >= fCperMIP_.size()) {
38  throw cms::Exception("OutOfBound") << "Trying to access thickness index " << thickness
39  << " in fCperMIP, which is of size " << fCperMIP_.size();
40  }
41  if (thickness >= chargeCollectionEfficiency_.size()) {
42  throw cms::Exception("OutOfBound") << "Trying to access thickness index " << thickness
43  << " in chargeCollectionEfficiency, which is of size "
45  }
46 
47  /* get the hardware pT in ADC counts: */
48  int hwPt = trgCell.hwPt();
49 
50  // Convert ADC to charge in fC (in EE+FH) or in MIPs (in BH)
51  double amplitude = hwPt * lsb_;
52 
55  }
56 
57  /* convert the charge amplitude in MIP: */
58  double trgCellMipP = amplitude;
59 
60  if (fCperMIP_[thickness] > 0) {
61  trgCellMipP /= fCperMIP_[thickness];
62  }
63 
64  /* compute the transverse-mip */
65  double trgCellMipPt = trgCellMipP / std::cosh(trgCell.eta());
66 
67  /* setting pT [mip] */
68  trgCell.setMipPt(trgCellMipPt);
69 }
bool isSilicon(const DetId &) const
int thicknessIndex(const DetId &) const
std::vector< double > chargeCollectionEfficiency_
uint32_t detId() const
void setMipPt(double value)
Definition: DetId.h:17
int hwPt() const
Definition: L1Candidate.h:35
double eta() const final
momentum pseudorapidity

◆ calibrateMipTinGeV()

void HGCalTriggerCellCalibration::calibrateMipTinGeV ( l1t::HGCalTriggerCell trgCell) const

Definition at line 71 of file HGCalTriggerCellCalibration.cc.

References ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), dEdX_weights_, l1t::HGCalTriggerCell::detId(), reco::LeafCandidate::eta(), Exception, HGCalTriggerTools::isSilicon(), HGCalTriggerTools::layerWithOffset(), l1t::HGCalTriggerCell::mipPt(), reco::LeafCandidate::phi(), reco::LeafCandidate::setP4(), Calorimetry_cff::thickness, thicknessCorrection_, HGCalTriggerTools::thicknessIndex(), and triggerTools_.

Referenced by calibrateInGeV().

71  {
72  constexpr double MevToGeV(0.001);
73  double trgCellEt(0.);
74 
75  DetId trgdetid(trgCell.detId());
76  unsigned trgCellLayer = triggerTools_.layerWithOffset(trgdetid);
77  bool isSilicon = triggerTools_.isSilicon(trgdetid);
78  constexpr int kScintillatorIndex = 0;
79  unsigned thickness = isSilicon ? triggerTools_.thicknessIndex(trgdetid) : kScintillatorIndex;
80  if (thickness >= thicknessCorrection_.size()) {
81  throw cms::Exception("OutOfBound") << "Trying to access thickness index " << thickness
82  << " in thicknessCorrection, which is of size " << thicknessCorrection_.size();
83  }
84 
85  /* weight the amplitude by the absorber coefficient in MeV/mip + bring it in
86  * GeV */
87  trgCellEt = trgCell.mipPt() * MevToGeV;
88  trgCellEt *= dEdX_weights_.at(trgCellLayer);
89 
90  /* correct for the cell-thickness */
92  trgCellEt /= thicknessCorrection_[thickness];
93  }
94 
95  math::PtEtaPhiMLorentzVector calibP4(trgCellEt, trgCell.eta(), trgCell.phi(), 0.);
96  trgCell.setP4(calibP4);
97 }
bool isSilicon(const DetId &) const
int thicknessIndex(const DetId &) const
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
uint32_t detId() const
unsigned layerWithOffset(const DetId &) const
Definition: DetId.h:17
std::vector< double > thicknessCorrection_
double phi() const final
momentum azimuthal angle
void setP4(const LorentzVector &p4) final
set 4-momentum
double eta() const final
momentum pseudorapidity

◆ setGeometry()

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

Member Data Documentation

◆ chargeCollectionEfficiency_

std::vector<double> HGCalTriggerCellCalibration::chargeCollectionEfficiency_
private

Definition at line 20 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateInMipT(), and HGCalTriggerCellCalibration().

◆ dEdX_weights_

std::vector<double> HGCalTriggerCellCalibration::dEdX_weights_
private

Definition at line 22 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateMipTinGeV().

◆ fCperMIP_

std::vector<double> HGCalTriggerCellCalibration::fCperMIP_
private

Definition at line 19 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateInMipT(), and HGCalTriggerCellCalibration().

◆ lsb_

double HGCalTriggerCellCalibration::lsb_
private

Definition at line 18 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateInMipT().

◆ thicknessCorrection_

std::vector<double> HGCalTriggerCellCalibration::thicknessCorrection_
private

◆ triggerTools_

HGCalTriggerTools HGCalTriggerCellCalibration::triggerTools_
private

Definition at line 24 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateInMipT(), calibrateMipTinGeV(), and setGeometry().