CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCalTriggerCellCalibration.cc
Go to the documentation of this file.
2 
3 #include <cmath>
4 
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 }
31 
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 
53  if (chargeCollectionEfficiency_[thickness] > 0) {
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 }
70 
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 */
91  if (thicknessCorrection_[thickness] > 0) {
92  trgCellEt /= thicknessCorrection_[thickness];
93  }
94 
95  math::PtEtaPhiMLorentzVector calibP4(trgCellEt, trgCell.eta(), trgCell.phi(), 0.);
96  trgCell.setP4(calibP4);
97 }
98 
100  /* calibrate from ADC count to transverse mip */
101  calibrateInMipT(trgCell);
102 
103  /* calibrate from mip count to GeV */
104  calibrateMipTinGeV(trgCell);
105 }
int thicknessIndex(const DetId &) const
unsigned layerWithOffset(const DetId &) const
std::vector< double > chargeCollectionEfficiency_
double mipPt() const
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
void calibrateMipTinGeV(l1t::HGCalTriggerCell &) const
uint32_t detId() const
void setMipPt(double value)
Definition: DetId.h:17
bool isSilicon(const DetId &) const
void calibrateInGeV(l1t::HGCalTriggerCell &) const
int hwPt() const
Definition: L1Candidate.h:35
std::vector< double > thicknessCorrection_
Log< level::Warning, false > LogWarning
void calibrateInMipT(l1t::HGCalTriggerCell &) const
double phi() const final
momentum azimuthal angle
void setP4(const LorentzVector &p4) final
set 4-momentum
HGCalTriggerCellCalibration(const edm::ParameterSet &conf)
double eta() const final
momentum pseudorapidity