CMS 3D CMS Logo

HGCalTriggerCellCalibration.cc
Go to the documentation of this file.
2 
3 //class constructor
5 
6  LSB_ = beCodecConfig.getParameter<double>("cellLSB");
7  fCperMIP_ee_ = beCodecConfig.getParameter<std::vector<double>>("fCperMIPee");
8  fCperMIP_fh_ = beCodecConfig.getParameter<std::vector<double>>("fCperMIPfh");
9  dEdX_weights_ = beCodecConfig.getParameter<std::vector<double>>("dEdXweights");
10  thickCorr_ = beCodecConfig.getParameter<std::vector<double>>("thickCorr");
11 }
12 
13 
15 {
16 
17  HGCalDetId trgdetid( trgCell.detId() );
18  int subdet = trgdetid.subdetId();
19 
20  /* get the hardware pT in ADC counts: */
21  int hwPt = trgCell.hwPt();
22 
23  /* set the lowest signal bit and convert in charge amplitude: */
24  double amplitude = hwPt * LSB_;
25 
26  /* convert the charge amplitude in MIP: */
27  if( subdet == HGCEE ){
28  amplitude = amplitude / fCperMIP_ee_.at(cellThickness-1);
29  }else if( subdet == HGCHEF ){
30  amplitude = amplitude / fCperMIP_fh_.at(cellThickness-1);
31  }else if( subdet == HGCHEB ){
32  edm::LogWarning("DataNotFound") << "WARNING: the BH trgCells are not yet implemented";
33  }
34 
35  /* correct the charge amplitude for the sensor thickness */
36  double trgCellMipP = amplitude * thickCorr_.at( cellThickness-1 );
37  double trgCellMipPt = trgCellMipP/cosh( trgCell.eta() );
38 
39  /* setting pT [mip] */
40  trgCell.setMipPt( trgCellMipPt ) ;
41 }
42 
43 
45 {
46  const double MevToGeV(0.001);
47 
48  HGCalDetId trgdetid( trgCell.detId() );
49  int trgCellLayer = trgdetid.layer();
50  int subdet = trgdetid.subdetId();
51 
52  /* get the transverse momentum in mip units */
53  double mipP = trgCell.mipPt() * cosh( trgCell.eta() );
54 
55  if( subdet == HGCHEF ){
56  trgCellLayer = trgCellLayer + 28;
57  }
58 
59  //weight the amplitude by the absorber coefficient in MeV/mip + bring it in GeV
60  double trgCellE = mipP * dEdX_weights_.at(trgCellLayer) * MevToGeV;
61 
62  //assign the new energy to the four-vector of the trigger cell
63  math::PtEtaPhiMLorentzVector calibP4(trgCellE/cosh( trgCell.eta() ),
64  trgCell.eta(),
65  trgCell.phi(),
66  trgCell.p4().M() );
67 
68  // overwriting the 4p with the calibrated 4p
69  trgCell.setP4( calibP4 );
70 
71 }
72 
74 {
75 
76  /* calibrate from ADC count to transverse mip */
77  calibrateInMipT(trgCell, cellThickness);
78 
79  /* calibrate from mip count to GeV */
80  calibrateMipTinGeV(trgCell);
81 
82 }
83 
T getParameter(std::string const &) const
virtual double eta() const final
momentum pseudorapidity
void calibrateInMipT(l1t::HGCalTriggerCell &, int cellThickness)
void calibrateInGeV(l1t::HGCalTriggerCell &, int cellThickness)
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 calibrateMipTinGeV(l1t::HGCalTriggerCell &)
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
int hwPt() const
Definition: L1Candidate.h:48
virtual void setP4(const LorentzVector &p4) final
set 4-momentum
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
int layer() const
get the layer #
Definition: HGCalDetId.h:48
HGCalTriggerCellCalibration(const edm::ParameterSet &conf)