CMS 3D CMS Logo

EgammaHLTHcalVarProducerFromRecHit.cc
Go to the documentation of this file.
1 // Class: EgammaHLTHcalVarProducerFromRecHit
2 
3 /*
4 
5 Author: Swagata Mukherjee
6 
7 Date: August 2021
8 
9 This class is similar to the existing class EgammaHLTBcHcalIsolationProducersRegional,
10 but the new feature in this code is that the HCAL recHits are used instead of the
11 calotowers which is expected to be phased out sometime in Run3.
12 The old class can also be used until calotowers stay. After that, one need to switch to this new one.
13 
14 As the old producer code, this one also produces either Hcal isolation or H for H/E depending if doEtSum=true or false.
15 H for H/E = either a single HCAL tower behind SC, or towers in a cone, and hcal isolation has these tower(s) excluded.
16 A rho correction can be applied
17 
18 */
19 
38 
40 public:
42 
43 public:
44  void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const final;
45  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
46 
47 private:
48  const bool doEtSum_;
53  const double innerCone_;
54  const double outerCone_;
55  const int depth_;
56  const int maxSeverityHB_;
57  const int maxSeverityHE_;
58  const bool useSingleTower_;
59  const bool doRhoCorrection_;
60  const double rhoScale_;
61  const double rhoMax_;
62  const std::vector<double> effectiveAreas_;
63  const std::vector<double> absEtaLowEdges_;
73 
74  //Get HCAL thresholds from GT
76  bool cutsFromDB;
77 };
78 
80  : doEtSum_(config.getParameter<bool>("doEtSum")),
81  eThresHB_(config.getParameter<EgammaHcalIsolation::arrayHB>("eThresHB")),
82  etThresHB_(config.getParameter<EgammaHcalIsolation::arrayHB>("etThresHB")),
83  eThresHE_(config.getParameter<EgammaHcalIsolation::arrayHE>("eThresHE")),
84  etThresHE_(config.getParameter<EgammaHcalIsolation::arrayHE>("etThresHE")),
85  innerCone_(config.getParameter<double>("innerCone")),
86  outerCone_(config.getParameter<double>("outerCone")),
87  depth_(config.getParameter<int>("depth")),
88  maxSeverityHB_(config.getParameter<int>("maxSeverityHB")),
89  maxSeverityHE_(config.getParameter<int>("maxSeverityHE")),
90  useSingleTower_(config.getParameter<bool>("useSingleTower")),
91  doRhoCorrection_(config.getParameter<bool>("doRhoCorrection")),
92  rhoScale_(config.getParameter<double>("rhoScale")),
93  rhoMax_(config.getParameter<double>("rhoMax")),
94  effectiveAreas_(config.getParameter<std::vector<double> >("effectiveAreas")),
95  absEtaLowEdges_(config.getParameter<std::vector<double> >("absEtaLowEdges")),
96  recoEcalCandidateProducer_(consumes(config.getParameter<edm::InputTag>("recoEcalCandidateProducer"))),
97  hbheRecHitsTag_(consumes(config.getParameter<edm::InputTag>("hbheRecHitsTag"))),
98  rhoProducer_(doRhoCorrection_ ? consumes<double>(config.getParameter<edm::InputTag>("rhoProducer"))
99  : edm::EDGetTokenT<double>()),
100  caloGeometryToken_{esConsumes()},
101  hcalTopologyToken_{esConsumes()},
102  hcalChannelQualityToken_{esConsumes(edm::ESInputTag("", "withTopo"))},
103  hcalSevLvlComputerToken_{esConsumes()},
104  caloTowerConstituentsMapToken_{esConsumes()},
105  putToken_{produces<reco::RecoEcalCandidateIsolationMap>()},
106  cutsFromDB(
107  config.getParameter<bool>("usePFThresholdsFromDB")) { //Retrieve HCAL PF thresholds - from config or from DB
108  if (doRhoCorrection_) {
109  if (absEtaLowEdges_.size() != effectiveAreas_.size()) {
110  throw cms::Exception("IncompatibleVects") << "absEtaLowEdges and effectiveAreas should be of the same size. \n";
111  }
112 
113  if (absEtaLowEdges_.at(0) != 0.0) {
114  throw cms::Exception("IncompleteCoverage") << "absEtaLowEdges should start from 0. \n";
115  }
116 
117  for (unsigned int aIt = 0; aIt < absEtaLowEdges_.size() - 1; aIt++) {
118  if (!(absEtaLowEdges_.at(aIt) < absEtaLowEdges_.at(aIt + 1))) {
119  throw cms::Exception("ImproperBinning") << "absEtaLowEdges entries should be in increasing order. \n";
120  }
121  }
122  }
123 
124  if (cutsFromDB) {
125  hcalCutsToken_ = esConsumes<HcalPFCuts, HcalPFCutsRcd>(edm::ESInputTag("", "withTopo"));
126  }
127 }
128 
131 
132  desc.add<edm::InputTag>("recoEcalCandidateProducer", edm::InputTag("hltRecoEcalCandidate"));
133  desc.add<edm::InputTag>("rhoProducer", edm::InputTag("fixedGridRhoFastjetAllCalo"));
134  desc.add<edm::InputTag>("hbheRecHitsTag", edm::InputTag("hltHbhereco"));
135  desc.add<bool>("doRhoCorrection", false);
136  desc.add<double>("rhoMax", 999999.);
137  desc.add<double>(("rhoScale"), 1.0);
138  //eThresHB/HE are from RecoParticleFlow/PFClusterProducer/python/particleFlowRecHitHBHE_cfi.py
139  desc.add<std::vector<double> >("eThresHB", {0.1, 0.2, 0.3, 0.3});
140  desc.add<std::vector<double> >("etThresHB", {0, 0, 0, 0});
141  desc.add<std::vector<double> >("eThresHE", {0.1, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2});
142  desc.add<std::vector<double> >("etThresHE", {0, 0, 0, 0, 0, 0, 0});
143  desc.add<bool>("usePFThresholdsFromDB", true);
144  desc.add<double>("innerCone", 0);
145  desc.add<double>("outerCone", 0.14);
146  desc.add<int>("depth", 0);
147  desc.add<int>("maxSeverityHB", 9);
148  desc.add<int>("maxSeverityHE", 9);
149  desc.add<bool>("doEtSum", false);
150  desc.add<bool>("useSingleTower", false);
151  desc.add<std::vector<double> >("effectiveAreas", {0.079, 0.25}); // 2016 post-ichep sinEle default
152  desc.add<std::vector<double> >("absEtaLowEdges", {0.0, 1.479}); // Barrel, Endcap
153  descriptions.add("hltEgammaHLTHcalVarProducerFromRecHit", desc);
154 }
155 
158  const edm::EventSetup &iSetup) const {
159  auto recoEcalCandHandle = iEvent.getHandle(recoEcalCandidateProducer_);
160 
161  double rho = 0.0;
162 
163  if (doRhoCorrection_) {
164  rho = iEvent.get(rhoProducer_);
165  if (rho > rhoMax_) {
166  rho = rhoMax_;
167  }
168  rho = rho * rhoScale_;
169  }
170 
171  reco::RecoEcalCandidateIsolationMap isoMap(recoEcalCandHandle);
172 
173  for (unsigned int iRecoEcalCand = 0; iRecoEcalCand < recoEcalCandHandle->size(); iRecoEcalCand++) {
174  reco::RecoEcalCandidateRef recoEcalCandRef(recoEcalCandHandle, iRecoEcalCand);
175 
176  float isol = 0;
179 
180  if (useSingleTower_) {
181  if (!doEtSum_) { //this is single tower based H/E
184  } else { //this is cone-based HCAL isolation with single tower based footprint removal
187  }
188  } else { //useSingleTower_=False means H/E is cone-based
191  }
192 
194  outerCone_,
195  internal,
196  innerCone_,
197  eThresHB_,
198  etThresHB_,
200  eThresHE_,
201  etThresHE_,
203  iEvent.get(hbheRecHitsTag_),
204  iSetup.getData(caloGeometryToken_),
205  iSetup.getData(hcalTopologyToken_),
209  const HcalPFCuts *hcalCuts{nullptr};
210  if (cutsFromDB) {
211  hcalCuts = &iSetup.getData(hcalCutsToken_);
212  }
213 
214  if (useSingleTower_) {
215  if (doEtSum_) { //this is cone-based HCAL isolation with single tower based footprint removal
216  isol = thisHcalVar_.getHcalEtSumBc(recoEcalCandRef.get(), depth_, hcalCuts); //depth=0 means all depths
217  } else { //this is single tower based H/E
218  isol = thisHcalVar_.getHcalESumBc(recoEcalCandRef.get(), depth_, hcalCuts); //depth=0 means all depths
219  }
220  } else { //useSingleTower_=False means H/E is cone-based.
221  if (doEtSum_) { //hcal iso
222  isol = thisHcalVar_.getHcalEtSum(recoEcalCandRef.get(), depth_, hcalCuts); //depth=0 means all depths
223  } else { // doEtSum_=False means sum up energy, this is for H/E
224  isol = thisHcalVar_.getHcalESum(recoEcalCandRef.get(), depth_, hcalCuts); //depth=0 means all depths
225  }
226  }
227 
228  if (doRhoCorrection_) {
229  int iEA = -1;
230  auto scEta = std::abs(recoEcalCandRef->superCluster()->eta());
231  for (int bIt = absEtaLowEdges_.size() - 1; bIt > -1; bIt--) {
232  if (scEta >= absEtaLowEdges_[bIt]) {
233  iEA = bIt;
234  break;
235  }
236  }
237  isol = isol - rho * effectiveAreas_[iEA];
238  }
239 
240  isoMap.insert(recoEcalCandRef, isol);
241  }
242 
243  iEvent.emplace(putToken_, isoMap);
244 }
245 
const edm::ESGetToken< HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd > hcalSevLvlComputerToken_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_
double getHcalESumBc(const reco::Candidate *c, int depth, const HcalPFCuts *hcalCuts) const
Definition: config.py:1
double getHcalEtSumBc(const reco::Candidate *c, int depth, const HcalPFCuts *hcalCuts) const
int iEvent
Definition: GenABIO.cc:224
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const final
const edm::ESGetToken< CaloTowerConstituentsMap, CaloGeometryRecord > caloTowerConstituentsMapToken_
const edm::ESGetToken< HcalChannelQuality, HcalChannelQualityRcd > hcalChannelQualityToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double getHcalESum(const reco::Candidate *c, int depth, const HcalPFCuts *hcalCuts) const
const edm::EDGetTokenT< HBHERecHitCollection > hbheRecHitsTag_
const edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
void insert(const key_type &k, const data_type &v)
insert an association
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::ESGetToken< HcalPFCuts, HcalPFCutsRcd > hcalCutsToken_
HLT enums.
double getHcalEtSum(const reco::Candidate *c, int depth, const HcalPFCuts *hcalCuts) const
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:226
std::array< double, 4 > arrayHB
const edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > hcalTopologyToken_
const edm::EDPutTokenT< reco::RecoEcalCandidateIsolationMap > putToken_
std::array< double, 7 > arrayHE