CMS 3D CMS Logo

HGCalTriggerCellCalibration.cc
Go to the documentation of this file.
2 
3 //class constructor
5  : LSB_silicon_fC_(conf.getParameter<double>("siliconCellLSB_fC")),
6  LSB_scintillator_MIP_(conf.getParameter<double>("scintillatorCellLSB_MIP")),
7  fCperMIP_(conf.getParameter<double>("fCperMIP")),
8  thickCorr_(conf.getParameter<double>("thickCorr")),
9  dEdX_weights_(conf.getParameter<std::vector<double>>("dEdXweights")) {
10  if (fCperMIP_ <= 0) {
11  edm::LogWarning("DivisionByZero") << "WARNING: the MIP->fC correction factor is zero or negative. It won't be "
12  "applied to correct trigger cell energies.";
13  }
14  if (thickCorr_ <= 0) {
15  edm::LogWarning("DivisionByZero") << "WARNING: the cell-thickness correction factor is zero or negative. It won't "
16  "be applied to correct trigger cell energies.";
17  }
18 }
19 
21  bool isSilicon = triggerTools_.isSilicon(trgCell.detId());
22 
23  /* get the hardware pT in ADC counts: */
24  int hwPt = trgCell.hwPt();
25 
26  // Convert ADC to charge in fC (in EE+FH) or in MIPs (in BH)
27  double amplitude = hwPt * (!isSilicon ? LSB_scintillator_MIP_ : LSB_silicon_fC_);
28 
29  // The responses of the different cell thicknesses have been equalized
30  // to the 200um response in the front-end. So there is only one global
31  // fCperMIP and thickCorr here
32  /* convert the charge amplitude in MIP: */
33  double trgCellMipP = amplitude;
34  if (isSilicon && fCperMIP_ > 0) {
35  trgCellMipP /= fCperMIP_;
36  }
37 
38  /* compute the transverse-mip */
39  double trgCellMipPt = trgCellMipP / cosh(trgCell.eta());
40 
41  /* setting pT [mip] */
42  trgCell.setMipPt(trgCellMipPt);
43 }
44 
46  const double MevToGeV(0.001);
47 
48  DetId trgdetid(trgCell.detId());
49  unsigned trgCellLayer = triggerTools_.layerWithOffset(trgdetid);
50 
51  if (dEdX_weights_.at(trgCellLayer) == 0.) {
52  throw cms::Exception("BadConfiguration")
53  << "Trigger cell energy forced to 0 by calibration coefficients.\n"
54  << "The configuration should be changed. "
55  << "Discarded layers should be defined in hgcalTriggerGeometryESProducer.TriggerGeometry.DisconnectedLayers "
56  "and not with calibration coefficients = 0\n";
57  }
58 
59  /* weight the amplitude by the absorber coefficient in MeV/mip + bring it in GeV */
60  double trgCellEt = trgCell.mipPt() * dEdX_weights_.at(trgCellLayer) * MevToGeV;
61 
62  /* correct for the cell-thickness */
63  if (triggerTools_.isSilicon(trgdetid) && thickCorr_ > 0) {
64  trgCellEt /= thickCorr_;
65  }
66  /* assign the new energy to the four-vector of the trigger cell */
67  math::PtEtaPhiMLorentzVector calibP4(trgCellEt, trgCell.eta(), trgCell.phi(), 0.);
68 
69  /* overwriting the 4p with the calibrated 4p */
70  trgCell.setP4(calibP4);
71 }
72 
74  /* calibrate from ADC count to transverse mip */
75  calibrateInMipT(trgCell);
76 
77  /* calibrate from mip count to GeV */
78  calibrateMipTinGeV(trgCell);
79 }
double eta() const final
momentum pseudorapidity
unsigned layerWithOffset(const DetId &) const
double mipPt() const
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
uint32_t detId() const
void setMipPt(double value)
void calibrateInMipT(l1t::HGCalTriggerCell &)
void calibrateMipTinGeV(l1t::HGCalTriggerCell &)
Definition: DetId.h:18
bool isSilicon(const DetId &) const
int hwPt() const
Definition: L1Candidate.h:48
void calibrateInGeV(l1t::HGCalTriggerCell &)
double phi() const final
momentum azimuthal angle
void setP4(const LorentzVector &p4) final
set 4-momentum
HGCalTriggerCellCalibration(const edm::ParameterSet &conf)