CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
FixedGridRhoProducerFastjetFromRecHit.cc
Go to the documentation of this file.
1 /*
2 
3 Author: Swagata Mukherjee
4 
5 Date: November 2021
6 
7 This producer computes rho from ECAL and HCAL recHits.
8 The current plan is to use it in egamma and muon HLT, who currently use
9 the other producer called FixedGridRhoProducerFastjet.
10 At HLT, FixedGridRhoProducerFastjet takes calotowers as input and computes rho.
11 But calotowers are expected to be phased out sometime in Run3.
12 So this recHit-based rho producer, FixedGridRhoProducerFastjetFromRecHit, can be used as an alternative.
13 
14 */
15 
28 #include "fastjet/tools/GridMedianBackgroundEstimator.hh"
30 
32 public:
35  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
36 
37 private:
38  void produce(edm::Event &, const edm::EventSetup &) override;
39  std::array<double, 4> getHitP4(const DetId &detId, const double hitE, const CaloGeometry &caloGeometry) const;
40  bool passedHcalNoiseCut(const HBHERecHit &hit) const;
41  bool passedEcalNoiseCut(const EcalRecHit &hit, const EcalPFRecHitThresholds &thresholds) const;
42 
43  fastjet::GridMedianBackgroundEstimator bge_;
47 
50 
51  //Muon HLT currently use ECAL-only rho for ECAL isolation,
52  //and HCAL-only rho for HCAL isolation. So, this skipping option is needed.
53  const bool skipHCAL_;
54  const bool skipECAL_;
55 
58 };
59 
61  : bge_(iConfig.getParameter<double>("maxRapidity"), iConfig.getParameter<double>("gridSpacing")),
62  hbheRecHitsTag_(consumes(iConfig.getParameter<edm::InputTag>("hbheRecHitsTag"))),
63  ebRecHitsTag_(consumes(iConfig.getParameter<edm::InputTag>("ebRecHitsTag"))),
64  eeRecHitsTag_(consumes(iConfig.getParameter<edm::InputTag>("eeRecHitsTag"))),
65  eThresHB_(iConfig.getParameter<EgammaHcalIsolation::arrayHB>("eThresHB")),
66  eThresHE_(iConfig.getParameter<EgammaHcalIsolation::arrayHE>("eThresHE")),
67  skipHCAL_(iConfig.getParameter<bool>("skipHCAL")),
68  skipECAL_(iConfig.getParameter<bool>("skipECAL")),
69  ecalPFRecHitThresholdsToken_{esConsumes()},
70  caloGeometryToken_{esConsumes()} {
71  if (skipHCAL_ && skipECAL_) {
72  throw cms::Exception("FixedGridRhoProducerFastjetFromRecHit")
73  << "skipHCAL and skipECAL both can't be True. Please make at least one of them False.";
74  }
75  produces<double>();
76 }
77 
80  //We use raw recHits and not the PF recHits, because in Phase1 muon and egamma HLT,
81  //PF recHit collection is regional, while the full raw recHit collection is available.
82  desc.add<edm::InputTag>("hbheRecHitsTag", edm::InputTag("hltHbhereco"));
83  desc.add<edm::InputTag>("ebRecHitsTag", edm::InputTag("hltEcalRecHit", "EcalRecHitsEB"));
84  desc.add<edm::InputTag>("eeRecHitsTag", edm::InputTag("hltEcalRecHit", "EcalRecHitsEE"));
85  desc.add<bool>("skipHCAL", false);
86  desc.add<bool>("skipECAL", false);
87  //eThresHB/HE are from RecoParticleFlow/PFClusterProducer/python/particleFlowRecHitHBHE_cfi.py
88  desc.add<std::vector<double> >("eThresHB", {0.1, 0.2, 0.3, 0.3});
89  desc.add<std::vector<double> >("eThresHE", {0.1, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2});
90  desc.add<double>("maxRapidity", 2.5);
91  desc.add<double>("gridSpacing", 0.55);
92  descriptions.addWithDefaultLabel(desc);
93 }
94 
96 
98  std::vector<fastjet::PseudoJet> inputs;
99  auto const &thresholds = iSetup.getData(ecalPFRecHitThresholdsToken_);
100  auto const &caloGeometry = iSetup.getData(caloGeometryToken_);
101 
102  if (!skipHCAL_) {
103  auto const &hbheRecHits = iEvent.get(hbheRecHitsTag_);
104  inputs.reserve(inputs.size() + hbheRecHits.size());
105  for (const auto &hit : hbheRecHits) {
106  if (passedHcalNoiseCut(hit)) {
107  const auto &hitp4 = getHitP4(hit.id(), hit.energy(), caloGeometry);
108  inputs.emplace_back(fastjet::PseudoJet(hitp4[0], hitp4[1], hitp4[2], hitp4[3]));
109  }
110  }
111  }
112 
113  if (!skipECAL_) {
114  auto const &ebRecHits = iEvent.get(ebRecHitsTag_);
115  inputs.reserve(inputs.size() + ebRecHits.size());
116  for (const auto &hit : ebRecHits) {
118  const auto &hitp4 = getHitP4(hit.id(), hit.energy(), caloGeometry);
119  inputs.emplace_back(fastjet::PseudoJet(hitp4[0], hitp4[1], hitp4[2], hitp4[3]));
120  }
121  }
122 
123  auto const &eeRecHits = iEvent.get(eeRecHitsTag_);
124  inputs.reserve(inputs.size() + eeRecHits.size());
125  for (const auto &hit : eeRecHits) {
127  const auto &hitp4 = getHitP4(hit.id(), hit.energy(), caloGeometry);
128  inputs.emplace_back(fastjet::PseudoJet(hitp4[0], hitp4[1], hitp4[2], hitp4[3]));
129  }
130  }
131  }
132 
133  bge_.set_particles(inputs);
134  iEvent.put(std::make_unique<double>(bge_.rho()));
135 }
136 
137 std::array<double, 4> FixedGridRhoProducerFastjetFromRecHit::getHitP4(const DetId &detId,
138  const double hitE,
139  const CaloGeometry &caloGeometry) const {
140  const CaloSubdetectorGeometry *subDetGeom = caloGeometry.getSubdetectorGeometry(detId);
141  const auto &gpPos = subDetGeom->getGeometry(detId)->repPos();
142  const double thispt = hitE / cosh(gpPos.eta());
143  const double thispx = thispt * cos(gpPos.phi());
144  const double thispy = thispt * sin(gpPos.phi());
145  const double thispz = thispt * sinh(gpPos.eta());
146  return std::array<double, 4>{{thispx, thispy, thispz, hitE}};
147 }
148 
149 //HCAL noise cleaning cuts.
151  const auto thisDetId = hit.id();
152  const auto thisDepth = thisDetId.depth();
153  return (thisDetId.subdet() == HcalBarrel && hit.energy() > eThresHB_[thisDepth - 1]) ||
154  (thisDetId.subdet() == HcalEndcap && hit.energy() > eThresHE_[thisDepth - 1]);
155 }
156 
157 //ECAL noise cleaning cuts using per-crystal PF-recHit thresholds.
159  const EcalPFRecHitThresholds &thresholds) const {
160  return (hit.energy() > thresholds[hit.detid()]);
161 }
162 
constexpr float energy() const
Definition: CaloRecHit.h:29
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const edm::EDGetTokenT< EcalRecHitCollection > ebRecHitsTag_
constexpr SubDetector subDetGeom[21]
std::array< double, 4 > getHitP4(const DetId &detId, const double hitE, const CaloGeometry &caloGeometry) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &, const edm::EventSetup &) override
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const DetId & detid() const
Definition: EcalRecHit.h:72
const edm::ESGetToken< EcalPFRecHitThresholds, EcalPFRecHitThresholdsRcd > ecalPFRecHitThresholdsToken_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
int iEvent
Definition: GenABIO.cc:224
constexpr HcalDetId id() const
get the id
Definition: HBHERecHit.h:41
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
float energy() const
Definition: EcalRecHit.h:68
ParameterDescriptionBase * add(U const &iLabel, T const &value)
unsigned int id
const edm::EDGetTokenT< HBHERecHitCollection > hbheRecHitsTag_
Definition: DetId.h:17
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
const edm::EDGetTokenT< EcalRecHitCollection > eeRecHitsTag_
FixedGridRhoProducerFastjetFromRecHit(const edm::ParameterSet &iConfig)
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_
bool passedEcalNoiseCut(const EcalRecHit &hit, const EcalPFRecHitThresholds &thresholds) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164
std::array< double, 4 > arrayHB
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::array< double, 7 > arrayHE