CMS 3D CMS Logo

HGCalRecHitWorkerSimple.cc
Go to the documentation of this file.
10 
13 {
15  tools_.reset(new hgcal::RecHitTools());
16  constexpr float keV2GeV = 1e-6;
17  // HGCee constants
18  hgcEE_keV2DIGI_ = ps.getParameter<double>("HGCEE_keV2DIGI");
19  hgcEE_fCPerMIP_ = ps.getParameter < std::vector<double> > ("HGCEE_fCPerMIP");
20  hgcEE_isSiFE_ = ps.getParameter<bool>("HGCEE_isSiFE");
22 
23  // HGChef constants
24  hgcHEF_keV2DIGI_ = ps.getParameter<double>("HGCHEF_keV2DIGI");
25  hgcHEF_fCPerMIP_ = ps.getParameter < std::vector<double> > ("HGCHEF_fCPerMIP");
26  hgcHEF_isSiFE_ = ps.getParameter<bool>("HGCHEF_isSiFE");
28 
29  // HGCheb constants
30  hgcHEB_keV2DIGI_ = ps.getParameter<double>("HGCHEB_keV2DIGI");
31  hgcHEB_isSiFE_ = ps.getParameter<bool>("HGCHEB_isSiFE");
33 
34  // layer weights (from Valeri/Arabella)
35  const auto& dweights = ps.getParameter < std::vector<double> > ("layerWeights");
36  for (auto weight : dweights)
37  {
38  weights_.push_back(weight);
39  }
40 
41  rechitMaker_->setLayerWeights(weights_);
42 
43  // residual correction for cell thickness
44  const auto& rcorr = ps.getParameter < std::vector<double> > ("thicknessCorrection");
45  rcorr_.clear();
46  rcorr_.push_back(1.f);
47  for (auto corr : rcorr)
48  {
49  rcorr_.push_back(1.0 / corr);
50  }
51 
52 
53  hgcEE_noise_fC_ = ps.getParameter < std::vector<double> > ("HGCEE_noise_fC");
54  hgcEE_cce_ = ps.getParameter< std::vector<double> > ("HGCEE_cce");
55  hgcHEF_noise_fC_ = ps.getParameter < std::vector<double> > ("HGCHEF_noise_fC");
56  hgcHEF_cce_ = ps.getParameter< std::vector<double> > ("HGCHEF_cce");
57  hgcHEB_noise_MIP_ = ps.getParameter<double>("HGCHEB_noise_MIP");
58 
59  // don't produce rechit if detid is a ghost one
60  rangeMatch_ = ps.getParameter<uint32_t>("rangeMatch");
61  rangeMask_ = ps.getParameter<uint32_t>("rangeMask");
62 }
63 
65 {
66  tools_->getEventSetup(es);
67  if (hgcEE_isSiFE_)
68  {
69  edm::ESHandle < HGCalGeometry > hgceeGeoHandle;
70  es.get<IdealGeometryRecord>().get("HGCalEESensitive", hgceeGeoHandle);
71  ddds_[0] = &(hgceeGeoHandle->topology().dddConstants());
72  }
73  else
74  {
75  ddds_[0] = nullptr;
76  }
77  if (hgcHEF_isSiFE_)
78  {
79  edm::ESHandle < HGCalGeometry > hgchefGeoHandle;
80  es.get<IdealGeometryRecord>().get("HGCalHESiliconSensitive", hgchefGeoHandle);
81  ddds_[1] = &(hgchefGeoHandle->topology().dddConstants());
82  }
83  else
84  {
85  ddds_[1] = nullptr;
86  }
87  ddds_[2] = nullptr;
88 }
89 
92 {
93  DetId detid = uncalibRH.id();
94 // don't produce rechit if detid is a ghost one
95 
96  if((detid & rangeMask_) == rangeMatch_)
97  return false;
98 
99  int thickness = -1;
100  float sigmaNoiseGeV = 0.f;
101  unsigned int layer = tools_->getLayerWithOffset(detid);
102  HGCalDetId hid;
103  float cce_correction = 1.0;
104 
105  switch (detid.subdetId())
106  {
107  case HGCEE:
108  rechitMaker_->setADCToGeVConstant(float(hgceeUncalib2GeV_));
109  hid = detid;
110  thickness = ddds_[hid.subdetId() - 3]->waferTypeL(hid.wafer());
111  cce_correction = hgcEE_cce_[thickness - 1];
112  sigmaNoiseGeV = 1e-3 * weights_[layer] * rcorr_[thickness]
113  * hgcEE_noise_fC_[thickness - 1] / hgcEE_fCPerMIP_[thickness - 1];
114  break;
115  case HGCHEF:
116  rechitMaker_->setADCToGeVConstant(float(hgchefUncalib2GeV_));
117  hid = detid;
118  thickness = ddds_[hid.subdetId() - 3]->waferTypeL(hid.wafer());
119  cce_correction = hgcHEF_cce_[thickness - 1];
120  sigmaNoiseGeV = 1e-3 * weights_[layer] * rcorr_[thickness]
121  * hgcHEF_noise_fC_[thickness - 1] / hgcHEF_fCPerMIP_[thickness - 1];
122  break;
123  case HcalEndcap:
124  case HGCHEB:
125  rechitMaker_->setADCToGeVConstant(float(hgchebUncalib2GeV_));
126  hid = detid;
127  sigmaNoiseGeV = 1e-3 * hgcHEB_noise_MIP_ * weights_[layer];
128  break;
129  default:
130  throw cms::Exception("NonHGCRecHit") << "Rechit with detid = " << detid.rawId()
131  << " is not HGC!";
132  }
133 
134  // make the rechit and put in the output collection
135 
136  HGCRecHit myrechit(rechitMaker_->makeRecHit(uncalibRH, 0));
137  const double new_E = myrechit.energy() * (thickness == -1 ? 1.0 : rcorr_[thickness])/cce_correction;
138 
139 
140  myrechit.setEnergy(new_E);
141  myrechit.setSignalOverSigmaNoise(new_E/sigmaNoiseGeV);
142  result.push_back(myrechit);
143 
144  return true;
145 }
146 
148 {
149 }
150 
T getParameter(std::string const &) const
std::vector< float > weights_
std::vector< double > hgcEE_cce_
std::array< const HGCalDDDConstants *, 3 > ddds_
std::vector< double > hgcHEF_fCPerMIP_
std::vector< double > hgcEE_fCPerMIP_
void push_back(T const &t)
Definition: weight.py:1
std::vector< double > hgcEE_noise_fC_
#define constexpr
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
std::vector< double > hgcHEF_cce_
float energy() const
Definition: CaloRecHit.h:17
std::vector< double > hgcHEF_noise_fC_
std::vector< double > rcorr_
const HGCalTopology & topology() const
Definition: HGCalGeometry.h:98
void set(const edm::EventSetup &es) override
double f[11][100]
int wafer() const
get the wafer #
Definition: HGCalDetId.h:42
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
JetCorrectorParameters corr
Definition: classes.h:5
std::unique_ptr< hgcal::RecHitTools > tools_
Definition: DetId.h:18
const HGCalDDDConstants & dddConstants() const
HGCalRecHitWorkerSimple(const edm::ParameterSet &)
const T & get() const
Definition: EventSetup.h:58
bool run(const edm::Event &evt, const HGCUncalibratedRecHit &uncalibRH, HGCRecHitCollection &result) override
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::unique_ptr< HGCalRecHitSimpleAlgo > rechitMaker_