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 &, int cellThickness)
 
void calibrateInMipT (l1t::HGCalTriggerCell &, int cellThickness)
 
void calibrateMipTinGeV (l1t::HGCalTriggerCell &)
 
 HGCalTriggerCellCalibration (const edm::ParameterSet &conf)
 
void print ()
 

Private Attributes

std::vector< double > dEdX_weights_
 
std::vector< double > fCperMIP_ee_
 
std::vector< double > fCperMIP_fh_
 
double LSB_
 
std::vector< double > thickCorr_
 

Detailed Description

Definition at line 17 of file HGCalTriggerCellCalibration.h.

Constructor & Destructor Documentation

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

Definition at line 4 of file HGCalTriggerCellCalibration.cc.

References corr, dEdX_weights_, fCperMIP_ee_, fCperMIP_fh_, edm::ParameterSet::getParameter(), LSB_, and thickCorr_.

4  {
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  for(auto corr : thickCorr_){
14  if(corr <= 0){
15  edm::LogWarning("DivisionByZero") << "WARNING: the cell-thickness correction factor is zero or negative. It won't be applied to correct trigger cell energies.";
16  }
17 
18  }
19 }
JetCorrectorParameters corr
Definition: classes.h:5

Member Function Documentation

void HGCalTriggerCellCalibration::calibrateInGeV ( l1t::HGCalTriggerCell trgCell,
int  cellThickness 
)

Definition at line 85 of file HGCalTriggerCellCalibration.cc.

References calibrateInMipT(), and calibrateMipTinGeV().

Referenced by HGCalTriggerBackend::HGCalTriggerSimCluster< FECODEC, DATA >::run().

86 {
87 
88  /* calibrate from ADC count to transverse mip */
89  calibrateInMipT(trgCell, cellThickness);
90 
91  /* calibrate from mip count to GeV */
92  calibrateMipTinGeV(trgCell);
93 
94 }
void calibrateInMipT(l1t::HGCalTriggerCell &, int cellThickness)
void calibrateMipTinGeV(l1t::HGCalTriggerCell &)
void HGCalTriggerCellCalibration::calibrateInMipT ( l1t::HGCalTriggerCell trgCell,
int  cellThickness 
)

Definition at line 22 of file HGCalTriggerCellCalibration.cc.

References CustomPhysics_cfi::amplitude, l1t::HGCalTriggerCell::detId(), reco::LeafCandidate::eta(), fCperMIP_ee_, fCperMIP_fh_, HGCEE, HGCHEB, HGCHEF, l1t::L1Candidate::hwPt(), LSB_, l1t::HGCalTriggerCell::setMipPt(), DetId::subdetId(), and thickCorr_.

Referenced by calibrateInGeV().

23 {
24 
25  HGCalDetId trgdetid( trgCell.detId() );
26  int subdet = trgdetid.subdetId();
27 
28  /* get the hardware pT in ADC counts: */
29  int hwPt = trgCell.hwPt();
30 
31  /* set the lowest signal bit and convert in charge amplitude: */
32  double amplitude = hwPt * LSB_;
33 
34  /* convert the charge amplitude in MIP: */
35  if( subdet == HGCEE ){
36  amplitude = amplitude / fCperMIP_ee_.at(cellThickness-1);
37  }else if( subdet == HGCHEF ){
38  amplitude = amplitude / fCperMIP_fh_.at(cellThickness-1);
39  }else if( subdet == HGCHEB ){
40  edm::LogWarning("DataNotFound") << "WARNING: the BH trgCells are not yet implemented";
41  }
42 
43  /* correct the charge amplitude for the sensor thickness */
44  double trgCellMipP = amplitude;
45  if( thickCorr_.at( cellThickness-1 ) > 0 ){
46  trgCellMipP = trgCellMipP / thickCorr_.at( cellThickness-1 );
47  }
48 
49  double trgCellMipPt = trgCellMipP/cosh( trgCell.eta() );
50 
51  /* setting pT [mip] */
52  trgCell.setMipPt( trgCellMipPt ) ;
53 }
virtual double eta() const final
momentum pseudorapidity
uint32_t detId() const
void setMipPt(double value)
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
void HGCalTriggerCellCalibration::calibrateMipTinGeV ( l1t::HGCalTriggerCell trgCell)

Definition at line 56 of file HGCalTriggerCellCalibration.cc.

References dEdX_weights_, l1t::HGCalTriggerCell::detId(), reco::LeafCandidate::eta(), HGCHEF, HGCalDetId::layer(), l1t::HGCalTriggerCell::mipPt(), reco::LeafCandidate::p4(), reco::LeafCandidate::phi(), and reco::LeafCandidate::setP4().

Referenced by calibrateInGeV().

57 {
58  const double MevToGeV(0.001);
59 
60  HGCalDetId trgdetid( trgCell.detId() );
61  int trgCellLayer = trgdetid.layer();
62  int subdet = trgdetid.subdetId();
63 
64  /* get the transverse momentum in mip units */
65  double mipP = trgCell.mipPt() * cosh( trgCell.eta() );
66 
67  if( subdet == HGCHEF ){
68  trgCellLayer = trgCellLayer + 28;
69  }
70 
71  //weight the amplitude by the absorber coefficient in MeV/mip + bring it in GeV
72  double trgCellE = mipP * dEdX_weights_.at(trgCellLayer) * MevToGeV;
73 
74  //assign the new energy to the four-vector of the trigger cell
75  math::PtEtaPhiMLorentzVector calibP4(trgCellE/cosh( trgCell.eta() ),
76  trgCell.eta(),
77  trgCell.phi(),
78  trgCell.p4().M() );
79 
80  // overwriting the 4p with the calibrated 4p
81  trgCell.setP4( calibP4 );
82 
83 }
virtual double eta() const final
momentum pseudorapidity
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
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
void HGCalTriggerCellCalibration::print ( )

Member Data Documentation

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

Definition at line 30 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateInMipT(), and HGCalTriggerCellCalibration().

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

Definition at line 31 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateInMipT(), and HGCalTriggerCellCalibration().

double HGCalTriggerCellCalibration::LSB_
private

Definition at line 29 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateInMipT(), and HGCalTriggerCellCalibration().

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

Definition at line 33 of file HGCalTriggerCellCalibration.h.

Referenced by calibrateInMipT(), and HGCalTriggerCellCalibration().