CMS 3D CMS Logo

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