Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

RecHit-based rho producer for HLT #36157

Merged
merged 10 commits into from
Nov 26, 2021
Merged
Show file tree
Hide file tree
Changes from 6 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions RecoJets/JetProducers/BuildFile.xml
Original file line number Diff line number Diff line change
Expand Up @@ -10,6 +10,7 @@
<use name="SimDataFormats/Vertex"/>
<use name="SimDataFormats/Track"/>
<use name="RecoJets/JetAlgorithms"/>
<use name="CondFormats/EcalObjects"/>
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I guess this should go to RecoJets/JetProducers/plugins/BuildFile.xml instead.

<use name="JetMETCorrections/JetCorrector"/>
<use name="DataFormats/CastorReco"/>
<use name="CommonTools/MVAUtils"/>
Expand Down
168 changes: 168 additions & 0 deletions RecoJets/JetProducers/plugins/FixedGridRhoProducerFastjetFromRecHit.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,168 @@
/*

Author: Swagata Mukherjee

Date: November 2021

This producer computes rho from ECAL and HCAL recHits.
The current plan is to use it in egamma and muon HLT, who currently use
the other producer called FixedGridRhoProducerFastjet.
At HLT, FixedGridRhoProducerFastjet takes calotowers as input and computes rho.
But calotowers are expected to be phased out sometime in Run3.
So this recHit-based rho producer, FixedGridRhoProducerFastjetFromRecHit, can be used as an alternative.

*/

#include "FWCore/Framework/interface/stream/EDProducer.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "DataFormats/Math/interface/LorentzVector.h"
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think this is no longer needed. Can you check the other headers as well, if there was anything left over that's no longer used?

#include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
#include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
#include "Geometry/Records/interface/CaloGeometryRecord.h"
#include "CondFormats/EcalObjects/interface/EcalPFRecHitThresholds.h"
#include "CondFormats/DataRecord/interface/EcalPFRecHitThresholdsRcd.h"
#include "fastjet/tools/GridMedianBackgroundEstimator.hh"
#include "RecoEgamma/EgammaIsolationAlgos/interface/EgammaHcalIsolation.h"

class FixedGridRhoProducerFastjetFromRecHit : public edm::stream::EDProducer<> {
public:
explicit FixedGridRhoProducerFastjetFromRecHit(const edm::ParameterSet &iConfig);
~FixedGridRhoProducerFastjetFromRecHit() override;
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);

private:
void produce(edm::Event &, const edm::EventSetup &) override;
std::array<double, 4> getHitP4(const DetId &detId, double hitE, const CaloGeometry &caloGeometry) const;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
std::array<double, 4> getHitP4(const DetId &detId, double hitE, const CaloGeometry &caloGeometry) const;
std::array<double, 4> getHitP4(const DetId &detId, const double hitE, const CaloGeometry &caloGeometry) const;

bool passedHcalNoiseCut(const HBHERecHit &hit) const;
bool passedEcalNoiseCut(const EcalRecHit &hit, const EcalPFRecHitThresholds &thresholds) const;

fastjet::GridMedianBackgroundEstimator bge_;
const edm::EDGetTokenT<HBHERecHitCollection> hbheRecHitsTag_;
const edm::EDGetTokenT<EcalRecHitCollection> ebRecHitsTag_;
const edm::EDGetTokenT<EcalRecHitCollection> eeRecHitsTag_;

const EgammaHcalIsolation::arrayHB eThresHB_;
const EgammaHcalIsolation::arrayHE eThresHE_;

//Muon HLT currently use ECAL-only rho for ECAL isolation,
//and HCAL-only rho for HCAL isolation. So, this skipping option is needed.
bool skipHCAL_;
bool skipECAL_;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
bool skipHCAL_;
bool skipECAL_;
const bool skipHCAL_;
const bool skipECAL_;


const edm::ESGetToken<EcalPFRecHitThresholds, EcalPFRecHitThresholdsRcd> ecalPFRecHitThresholdsToken_;
const edm::ESGetToken<CaloGeometry, CaloGeometryRecord> caloGeometryToken_;
};

FixedGridRhoProducerFastjetFromRecHit::FixedGridRhoProducerFastjetFromRecHit(const edm::ParameterSet &iConfig)
: bge_(iConfig.getParameter<double>("maxRapidity"), iConfig.getParameter<double>("gridSpacing")),
hbheRecHitsTag_(consumes(iConfig.getParameter<edm::InputTag>("hbheRecHitsTag"))),
ebRecHitsTag_(consumes(iConfig.getParameter<edm::InputTag>("ebRecHitsTag"))),
eeRecHitsTag_(consumes(iConfig.getParameter<edm::InputTag>("eeRecHitsTag"))),
eThresHB_(iConfig.getParameter<EgammaHcalIsolation::arrayHB>("eThresHB")),
eThresHE_(iConfig.getParameter<EgammaHcalIsolation::arrayHE>("eThresHE")),
Comment on lines +65 to +66
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I was wondering what happens if the input vector from the PSet does not have the size expected by arrayHB (4) and arrayHE (7). If this fails silently, it might be worth adding an explicitly check on the size.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

this is a good point, I am now checking what happens if I provide 2 elements instead for 4 for HB. It raise the following error and does not run; so I'm not adding any explicit check on size.

----- Begin Fatal Exception 18-Nov-2021 11:12:35 CET-----------------------
An exception of category 'Configuration' occurred while
   [0] Constructing the EventProcessor
   [1] Constructing module: class=FixedGridRhoProducerFastjetFromRecHit label='hltRhoNew'
Exception Message:
The parameter 'eThresHB' should have 4 elements, but has 2 elements in the configuration.
----- End Fatal Exception -------------------------------------------------

skipHCAL_(iConfig.getParameter<bool>("skipHCAL")),
skipECAL_(iConfig.getParameter<bool>("skipECAL")),
ecalPFRecHitThresholdsToken_{esConsumes()},
caloGeometryToken_{esConsumes()} {
if (skipHCAL_ && skipECAL_) {
throw cms::Exception("FixedGridRhoProducerFastjetFromRecHit")
<< "skipHCAL and skipECAL both can't be True. Please make at least one of them False.";
}
produces<double>();

This comment was marked as resolved.

}

void FixedGridRhoProducerFastjetFromRecHit::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
edm::ParameterSetDescription desc;
//We use raw recHits and not the PF recHits, because in Phase1 muon and egamma HLT,
//PF recHit collection is regional, while the full raw recHit collection is available.
desc.add<edm::InputTag>("hbheRecHitsTag", edm::InputTag("hltHbhereco"));
desc.add<edm::InputTag>("ebRecHitsTag", edm::InputTag("hltEcalRecHit", "EcalRecHitsEB"));
desc.add<edm::InputTag>("eeRecHitsTag", edm::InputTag("hltEcalRecHit", "EcalRecHitsEE"));
desc.add<bool>("skipHCAL", false);
desc.add<bool>("skipECAL", false);
//eThresHB/HE are from RecoParticleFlow/PFClusterProducer/python/particleFlowRecHitHBHE_cfi.py
desc.add<std::vector<double> >("eThresHB", {0.1, 0.2, 0.3, 0.3});
desc.add<std::vector<double> >("eThresHE", {0.1, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2});
desc.add<double>("maxRapidity", 2.5);
desc.add<double>("gridSpacing", 0.55);
descriptions.addWithDefaultLabel(desc);
}

FixedGridRhoProducerFastjetFromRecHit::~FixedGridRhoProducerFastjetFromRecHit() = default;

void FixedGridRhoProducerFastjetFromRecHit::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
std::vector<fastjet::PseudoJet> inputs;
auto const &thresholds = iSetup.getData(ecalPFRecHitThresholdsToken_);
auto const &caloGeometry = iSetup.getData(caloGeometryToken_);

This comment was marked as resolved.

if (!skipHCAL_) {
auto const &hbheRecHits = iEvent.get(hbheRecHitsTag_);
inputs.reserve(inputs.size() + hbheRecHits.size());
for (const auto &hit : hbheRecHits) {
if (passedHcalNoiseCut(hit)) {
const auto &hitp4 = getHitP4(hit.id(), hit.energy(), caloGeometry);
inputs.emplace_back(fastjet::PseudoJet(hitp4[0], hitp4[1], hitp4[2], hitp4[3]));
}
}
}

if (!skipECAL_) {
auto const &ebRecHits = iEvent.get(ebRecHitsTag_);
inputs.reserve(inputs.size() + ebRecHits.size());
for (const auto &hit : ebRecHits) {
if (passedEcalNoiseCut(hit, thresholds)) {
const auto &hitp4 = getHitP4(hit.id(), hit.energy(), caloGeometry);
inputs.emplace_back(fastjet::PseudoJet(hitp4[0], hitp4[1], hitp4[2], hitp4[3]));
}
}

auto const &eeRecHits = iEvent.get(eeRecHitsTag_);
inputs.reserve(inputs.size() + eeRecHits.size());
for (const auto &hit : eeRecHits) {
if (passedEcalNoiseCut(hit, thresholds)) {
const auto &hitp4 = getHitP4(hit.id(), hit.energy(), caloGeometry);
inputs.emplace_back(fastjet::PseudoJet(hitp4[0], hitp4[1], hitp4[2], hitp4[3]));
}
}
}

bge_.set_particles(inputs);
iEvent.put(std::make_unique<double>(bge_.rho()));
}

std::array<double, 4> FixedGridRhoProducerFastjetFromRecHit::getHitP4(const DetId &detId,
double hitE,
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
double hitE,
const double hitE,

const CaloGeometry &caloGeometry) const {
const CaloSubdetectorGeometry *subDetGeom = caloGeometry.getSubdetectorGeometry(detId);
const auto &gpPos = subDetGeom->getGeometry(detId)->repPos();
Comment on lines +140 to +141
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I take it that these pointers are always guaranteed to be valid (at least, the way they are used in this plugin). Is that correct?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

yes I think so

const double thispt = hitE / cosh(gpPos.eta());
const double thispx = thispt * cos(gpPos.phi());
const double thispy = thispt * sin(gpPos.phi());
const double thispz = thispt * sinh(gpPos.eta());
std::array<double, 4> hitp4{{thispx, thispy, thispz, hitE}};
return hitp4;
Copy link
Contributor

@VinInn VinInn Nov 21, 2021

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
return hitp4;
return std::array<double, 4> {{thispx, thispy, thispz, hitE}};

}

//HCAL noise cleaning cuts.
bool FixedGridRhoProducerFastjetFromRecHit::passedHcalNoiseCut(const HBHERecHit &hit) const {
const auto thisDetId = hit.id();
const auto thisDepth = thisDetId.depth();
if (thisDetId.subdet() == HcalBarrel && hit.energy() > eThresHB_[thisDepth - 1])
return true;
else if (thisDetId.subdet() == HcalEndcap && hit.energy() > eThresHE_[thisDepth - 1])
return true;
return false;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
return false;
return (thisDetId.subdet() == HcalBarrel && hit.energy() > eThresHB_[thisDepth - 1]) ||
(thisDetId.subdet() == HcalEndcap && hit.energy() > eThresHE_[thisDepth - 1]);

}

//ECAL noise cleaning cuts using per-crystal PF-recHit thresholds.
bool FixedGridRhoProducerFastjetFromRecHit::passedEcalNoiseCut(const EcalRecHit &hit,
const EcalPFRecHitThresholds &thresholds) const {
return (hit.energy() > thresholds[hit.detid()]);
}

DEFINE_FWK_MODULE(FixedGridRhoProducerFastjetFromRecHit);