CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Types | Protected Attributes
HGCalRecHitWorkerSimple Class Reference

#include <HGCalRecHitWorkerSimple.h>

Inheritance diagram for HGCalRecHitWorkerSimple:
HGCalRecHitWorkerBaseClass

Public Member Functions

 HGCalRecHitWorkerSimple (const edm::ParameterSet &)
 
bool run (const edm::Event &evt, const HGCUncalibratedRecHit &uncalibRH, HGCRecHitCollection &result) override
 
void set (const edm::EventSetup &es) override
 
 ~HGCalRecHitWorkerSimple () override
 
- Public Member Functions inherited from HGCalRecHitWorkerBaseClass
 HGCalRecHitWorkerBaseClass (const edm::ParameterSet &)
 
virtual ~HGCalRecHitWorkerBaseClass ()
 

Protected Types

enum  detectortype { hgcee =1, hgcfh =2, hgcbh =3 }
 

Protected Attributes

std::array< const HGCalDDDConstants *, 3 > ddds_
 
std::vector< double > hgcEE_cce_
 
std::vector< double > hgcEE_fCPerMIP_
 
bool hgcEE_isSiFE_
 
double hgcEE_keV2DIGI_
 
std::vector< double > hgcEE_noise_fC_
 
double hgceeUncalib2GeV_
 
bool hgcHEB_isSiFE_
 
double hgcHEB_keV2DIGI_
 
double hgcHEB_noise_MIP_
 
double hgchebUncalib2GeV_
 
std::vector< double > hgcHEF_cce_
 
std::vector< double > hgcHEF_fCPerMIP_
 
bool hgcHEF_isSiFE_
 
double hgcHEF_keV2DIGI_
 
std::vector< double > hgcHEF_noise_fC_
 
double hgchefUncalib2GeV_
 
bool killDeadChannels_
 
uint32_t rangeMask_
 
uint32_t rangeMatch_
 
std::vector< double > rcorr_
 
std::unique_ptr< HGCalRecHitSimpleAlgorechitMaker_
 
std::unique_ptr< hgcal::RecHitToolstools_
 
std::vector< int > v_chstatus_
 
std::vector< int > v_DB_reco_flags_
 
std::vector< float > weights_
 

Detailed Description

Definition at line 19 of file HGCalRecHitWorkerSimple.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

HGCalRecHitWorkerSimple::HGCalRecHitWorkerSimple ( const edm::ParameterSet ps)

Definition at line 9 of file HGCalRecHitWorkerSimple.cc.

References constexpr, corr, MillePedeFileConverter_cfg::e, f, edm::ParameterSet::getParameter(), hgcEE_cce_, hgcEE_fCPerMIP_, hgcEE_isSiFE_, hgcEE_keV2DIGI_, hgcEE_noise_fC_, hgceeUncalib2GeV_, hgcHEB_isSiFE_, hgcHEB_keV2DIGI_, hgcHEB_noise_MIP_, hgchebUncalib2GeV_, hgcHEF_cce_, hgcHEF_fCPerMIP_, hgcHEF_isSiFE_, hgcHEF_keV2DIGI_, hgcHEF_noise_fC_, hgchefUncalib2GeV_, rangeMask_, rangeMatch_, rcorr_, rechitMaker_, tools_, and weights_.

9  :
11 {
13  tools_.reset(new hgcal::RecHitTools());
14  constexpr float keV2GeV = 1e-6;
15  // HGCee constants
16  hgcEE_keV2DIGI_ = ps.getParameter<double>("HGCEE_keV2DIGI");
17  hgcEE_fCPerMIP_ = ps.getParameter < std::vector<double> > ("HGCEE_fCPerMIP");
18  hgcEE_isSiFE_ = ps.getParameter<bool>("HGCEE_isSiFE");
20 
21  // HGChef constants
22  hgcHEF_keV2DIGI_ = ps.getParameter<double>("HGCHEF_keV2DIGI");
23  hgcHEF_fCPerMIP_ = ps.getParameter < std::vector<double> > ("HGCHEF_fCPerMIP");
24  hgcHEF_isSiFE_ = ps.getParameter<bool>("HGCHEF_isSiFE");
26 
27  // HGCheb constants
28  hgcHEB_keV2DIGI_ = ps.getParameter<double>("HGCHEB_keV2DIGI");
29  hgcHEB_isSiFE_ = ps.getParameter<bool>("HGCHEB_isSiFE");
31 
32  // layer weights (from Valeri/Arabella)
33  const auto& dweights = ps.getParameter < std::vector<double> > ("layerWeights");
34  for (auto weight : dweights)
35  {
36  weights_.push_back(weight);
37  }
38 
39  rechitMaker_->setLayerWeights(weights_);
40 
41  // residual correction for cell thickness
42  const auto& rcorr = ps.getParameter < std::vector<double> > ("thicknessCorrection");
43  rcorr_.clear();
44  rcorr_.push_back(1.f);
45  for (auto corr : rcorr)
46  {
47  rcorr_.push_back(1.0 / corr);
48  }
49 
50 
51  hgcEE_noise_fC_ = ps.getParameter<edm::ParameterSet>("HGCEE_noise_fC").getParameter < std::vector<double> > ("values");
52  hgcEE_cce_ = ps.getParameter<edm::ParameterSet>("HGCEE_cce").getParameter< std::vector<double> > ("values");
53  hgcHEF_noise_fC_ = ps.getParameter<edm::ParameterSet>("HGCHEF_noise_fC").getParameter < std::vector<double> > ("values");
54  hgcHEF_cce_ = ps.getParameter<edm::ParameterSet>("HGCHEF_cce").getParameter< std::vector<double> > ("values");
55  hgcHEB_noise_MIP_ = ps.getParameter<edm::ParameterSet>("HGCHEB_noise_MIP").getParameter<double>("value");
56 
57  // don't produce rechit if detid is a ghost one
58  rangeMatch_ = ps.getParameter<uint32_t>("rangeMatch");
59  rangeMask_ = ps.getParameter<uint32_t>("rangeMask");
60 }
T getParameter(std::string const &) const
std::vector< float > weights_
std::vector< double > hgcEE_cce_
std::vector< double > hgcHEF_fCPerMIP_
std::vector< double > hgcEE_fCPerMIP_
Definition: weight.py:1
std::vector< double > hgcEE_noise_fC_
#define constexpr
std::vector< double > hgcHEF_cce_
std::vector< double > hgcHEF_noise_fC_
std::vector< double > rcorr_
double f[11][100]
HGCalRecHitWorkerBaseClass(const edm::ParameterSet &)
JetCorrectorParameters corr
Definition: classes.h:5
std::unique_ptr< hgcal::RecHitTools > tools_
std::unique_ptr< HGCalRecHitSimpleAlgo > rechitMaker_
HGCalRecHitWorkerSimple::~HGCalRecHitWorkerSimple ( )
override

Definition at line 169 of file HGCalRecHitWorkerSimple.cc.

References DEFINE_EDM_PLUGIN.

170 {
171 }

Member Function Documentation

bool HGCalRecHitWorkerSimple::run ( const edm::Event evt,
const HGCUncalibratedRecHit uncalibRH,
HGCRecHitCollection result 
)
overridevirtual

Implements HGCalRecHitWorkerBaseClass.

Definition at line 88 of file HGCalRecHitWorkerSimple.cc.

References ddds_, DetId::det(), MillePedeFileConverter_cfg::e, CaloRecHit::energy(), Exception, DetId::Forward, DetId::Hcal, HcalEndcap, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, hgcbh, HGCEE, hgcee, hgcEE_cce_, hgcEE_fCPerMIP_, hgcEE_noise_fC_, hgceeUncalib2GeV_, hgcfh, HGCHEB, hgcHEB_noise_MIP_, hgchebUncalib2GeV_, HGCHEF, hgcHEF_cce_, hgcHEF_fCPerMIP_, hgcHEF_noise_fC_, hgchefUncalib2GeV_, HGCUncalibratedRecHit::id(), edm::SortedCollection< T, SORT >::push_back(), rangeMask_, rangeMatch_, DetId::rawId(), rcorr_, rechitMaker_, DetId::subdetId(), tools_, HGCSiliconDetId::type(), and weights_.

90 {
91  DetId detid = uncalibRH.id();
92  // don't produce rechit if detid is a ghost one
93 
94  if (detid.det() == DetId::Forward || detid.det() == DetId::Hcal) {
95  if((detid & rangeMask_) == rangeMatch_) return false;
96  }
97 
98  int thickness = -1;
99  float sigmaNoiseGeV = 0.f;
100  unsigned int layer = tools_->getLayerWithOffset(detid);
101  float cce_correction = 1.0;
102  int idtype(0);
103 
104  switch (detid.det()) {
105  case DetId::HGCalEE:
106  idtype = hgcee;
107  thickness = 1 + HGCSiliconDetId(detid).type();
108  break;
109  case DetId::HGCalHSi:
110  idtype = hgcfh;
111  thickness = 1 + HGCSiliconDetId(detid).type();
112  break;
113  case DetId::HGCalHSc:
114  idtype = hgcbh;
115  break;
116  default:
117  switch (detid.subdetId()) {
118  case HGCEE:
119  idtype = hgcee;
120  thickness = ddds_[detid.subdetId()-3]->waferTypeL(HGCalDetId(detid).wafer());
121  break;
122  case HGCHEF:
123  idtype = hgcfh;
124  thickness = ddds_[detid.subdetId()-3]->waferTypeL(HGCalDetId(detid).wafer());
125  break;
126  case HcalEndcap:
127  case HGCHEB:
128  idtype = hgcbh;
129  break;
130  default:
131  break;
132  }
133  }
134  switch (idtype) {
135  case hgcee:
136  rechitMaker_->setADCToGeVConstant(float(hgceeUncalib2GeV_));
137  cce_correction = hgcEE_cce_[thickness - 1];
138  sigmaNoiseGeV = 1e-3 * weights_[layer] * rcorr_[thickness]
139  * hgcEE_noise_fC_[thickness - 1] / hgcEE_fCPerMIP_[thickness - 1];
140  break;
141  case hgcfh:
142  rechitMaker_->setADCToGeVConstant(float(hgchefUncalib2GeV_));
143  cce_correction = hgcHEF_cce_[thickness - 1];
144  sigmaNoiseGeV = 1e-3 * weights_[layer] * rcorr_[thickness]
145  * hgcHEF_noise_fC_[thickness - 1] / hgcHEF_fCPerMIP_[thickness - 1];
146  break;
147  case hgcbh:
148  rechitMaker_->setADCToGeVConstant(float(hgchebUncalib2GeV_));
149  sigmaNoiseGeV = 1e-3 * hgcHEB_noise_MIP_ * weights_[layer];
150  break;
151  default:
152  throw cms::Exception("NonHGCRecHit") << "Rechit with detid = "
153  << detid.rawId() << " is not HGC!";
154  }
155 
156  // make the rechit and put in the output collection
157 
158  HGCRecHit myrechit(rechitMaker_->makeRecHit(uncalibRH, 0));
159  const double new_E = myrechit.energy() * (thickness == -1 ? 1.0 : rcorr_[thickness])/cce_correction;
160 
161 
162  myrechit.setEnergy(new_E);
163  myrechit.setSignalOverSigmaNoise(new_E/sigmaNoiseGeV);
164  result.push_back(myrechit);
165 
166  return true;
167 }
constexpr float energy() const
Definition: CaloRecHit.h:31
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_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:47
void push_back(T const &t)
std::vector< double > hgcEE_noise_fC_
std::vector< double > hgcHEF_cce_
std::vector< double > hgcHEF_noise_fC_
int type() const
get the type
std::vector< double > rcorr_
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
std::unique_ptr< hgcal::RecHitTools > tools_
Definition: DetId.h:18
std::unique_ptr< HGCalRecHitSimpleAlgo > rechitMaker_
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39
void HGCalRecHitWorkerSimple::set ( const edm::EventSetup es)
overridevirtual

Implements HGCalRecHitWorkerBaseClass.

Definition at line 62 of file HGCalRecHitWorkerSimple.cc.

References HGCalTopology::dddConstants(), ddds_, edm::EventSetup::get(), hgcEE_isSiFE_, hgcHEF_isSiFE_, tools_, and HGCalGeometry::topology().

63 {
64  tools_->getEventSetup(es);
65  if (hgcEE_isSiFE_)
66  {
67  edm::ESHandle < HGCalGeometry > hgceeGeoHandle;
68  es.get<IdealGeometryRecord>().get("HGCalEESensitive", hgceeGeoHandle);
69  ddds_[0] = &(hgceeGeoHandle->topology().dddConstants());
70  }
71  else
72  {
73  ddds_[0] = nullptr;
74  }
75  if (hgcHEF_isSiFE_)
76  {
77  edm::ESHandle < HGCalGeometry > hgchefGeoHandle;
78  es.get<IdealGeometryRecord>().get("HGCalHESiliconSensitive", hgchefGeoHandle);
79  ddds_[1] = &(hgchefGeoHandle->topology().dddConstants());
80  }
81  else
82  {
83  ddds_[1] = nullptr;
84  }
85  ddds_[2] = nullptr;
86 }
std::array< const HGCalDDDConstants *, 3 > ddds_
const HGCalTopology & topology() const
std::unique_ptr< hgcal::RecHitTools > tools_
const HGCalDDDConstants & dddConstants() const
T get() const
Definition: EventSetup.h:63

Member Data Documentation

std::array<const HGCalDDDConstants*, 3> HGCalRecHitWorkerSimple::ddds_
protected

Definition at line 47 of file HGCalRecHitWorkerSimple.h.

Referenced by run(), and set().

std::vector<double> HGCalRecHitWorkerSimple::hgcEE_cce_
protected

Definition at line 33 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

std::vector<double> HGCalRecHitWorkerSimple::hgcEE_fCPerMIP_
protected

Definition at line 32 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

bool HGCalRecHitWorkerSimple::hgcEE_isSiFE_
protected

Definition at line 38 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and set().

double HGCalRecHitWorkerSimple::hgcEE_keV2DIGI_
protected

Definition at line 31 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple().

std::vector<double> HGCalRecHitWorkerSimple::hgcEE_noise_fC_
protected

Definition at line 42 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

double HGCalRecHitWorkerSimple::hgceeUncalib2GeV_
protected

Definition at line 31 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

bool HGCalRecHitWorkerSimple::hgcHEB_isSiFE_
protected

Definition at line 38 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple().

double HGCalRecHitWorkerSimple::hgcHEB_keV2DIGI_
protected

Definition at line 37 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple().

double HGCalRecHitWorkerSimple::hgcHEB_noise_MIP_
protected

Definition at line 44 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

double HGCalRecHitWorkerSimple::hgchebUncalib2GeV_
protected

Definition at line 37 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

std::vector<double> HGCalRecHitWorkerSimple::hgcHEF_cce_
protected

Definition at line 36 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

std::vector<double> HGCalRecHitWorkerSimple::hgcHEF_fCPerMIP_
protected

Definition at line 35 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

bool HGCalRecHitWorkerSimple::hgcHEF_isSiFE_
protected

Definition at line 38 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and set().

double HGCalRecHitWorkerSimple::hgcHEF_keV2DIGI_
protected

Definition at line 34 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple().

std::vector<double> HGCalRecHitWorkerSimple::hgcHEF_noise_fC_
protected

Definition at line 43 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

double HGCalRecHitWorkerSimple::hgchefUncalib2GeV_
protected

Definition at line 34 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

bool HGCalRecHitWorkerSimple::killDeadChannels_
protected

Definition at line 52 of file HGCalRecHitWorkerSimple.h.

uint32_t HGCalRecHitWorkerSimple::rangeMask_
protected

Definition at line 55 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

uint32_t HGCalRecHitWorkerSimple::rangeMatch_
protected

Definition at line 54 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

std::vector<double> HGCalRecHitWorkerSimple::rcorr_
protected

Definition at line 57 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

std::unique_ptr<HGCalRecHitSimpleAlgo> HGCalRecHitWorkerSimple::rechitMaker_
protected

Definition at line 59 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().

std::unique_ptr<hgcal::RecHitTools> HGCalRecHitWorkerSimple::tools_
protected

Definition at line 60 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), run(), and set().

std::vector<int> HGCalRecHitWorkerSimple::v_chstatus_
protected

Definition at line 49 of file HGCalRecHitWorkerSimple.h.

std::vector<int> HGCalRecHitWorkerSimple::v_DB_reco_flags_
protected

Definition at line 51 of file HGCalRecHitWorkerSimple.h.

std::vector<float> HGCalRecHitWorkerSimple::weights_
protected

Definition at line 58 of file HGCalRecHitWorkerSimple.h.

Referenced by HGCalRecHitWorkerSimple(), and run().