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

Simplify jet daughter processing in DeepBoostedJet #28035

Merged
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
20 changes: 4 additions & 16 deletions PhysicsTools/PatAlgos/python/tools/jetTools.py
Expand Up @@ -638,27 +638,15 @@ def setupBTagging(process, jetSource, pfCandidates, explicitJTA, pvSource, svSou
if btagInfo == 'pfDeepBoostedJetTagInfos':
if pfCandidates.value() == 'packedPFCandidates':
# case 1: running over jets whose daughters are PackedCandidates (only via updateJetCollection for now)
jetSrcName = jetSource.value().lower()
if 'updated' in jetSrcName:
puppi_value_map = ""
vertex_associator = ""
if 'withpuppidaughter' in jetSrcName:
# special case for Puppi jets reclustered from MiniAOD by analyzers
# need to specify 'WithPuppiDaughters' in the postfix when calling updateJetCollection
# daughters of these jets are already scaled by their puppi weights
has_puppi_weighted_daughters = True
else:
# default case for updating jet collection stored in MiniAOD, e.g., slimmedJetsAK8
# daughters are links to the original PackedCandidates, so NOT scaled by their puppi weights yet
has_puppi_weighted_daughters = False
else:
if 'updated' not in jetSource.value().lower():
raise ValueError("Invalid jet collection: %s. pfDeepBoostedJetTagInfos only supports running via updateJetCollection." % jetSource.value())
puppi_value_map = ""
vertex_associator = ""
elif pfCandidates.value() == 'particleFlow':
raise ValueError("Running pfDeepBoostedJetTagInfos with reco::PFCandidates is currently not supported.")
# case 2: running on new jet collection whose daughters are PFCandidates (e.g., cluster jets in RECO/AOD)
# daughters are the particles used in jet clustering, so already scaled by their puppi weights
# Uncomment the lines below after running pfDeepBoostedJetTagInfos with reco::PFCandidates becomes supported
# has_puppi_weighted_daughters = True
# puppi_value_map = "puppi"
# vertex_associator = "primaryVertexAssociation:original"
else:
Expand All @@ -668,7 +656,7 @@ def setupBTagging(process, jetSource, pfCandidates, explicitJTA, pvSource, svSou
jets = jetSource,
vertices = pvSource,
secondary_vertices = svSource,
has_puppi_weighted_daughters = has_puppi_weighted_daughters,
pf_candidates = pfCandidates,
puppi_value_map = puppi_value_map,
vertex_associator = vertex_associator,
),
Expand Down
38 changes: 16 additions & 22 deletions RecoBTag/FeatureTools/plugins/DeepBoostedJetTagInfoProducer.cc
Expand Up @@ -34,6 +34,7 @@ class DeepBoostedJetTagInfoProducer : public edm::stream::EDProducer<> {
typedef std::vector<reco::DeepBoostedJetTagInfo> DeepBoostedJetTagInfoCollection;
typedef reco::VertexCompositePtrCandidateCollection SVCollection;
typedef reco::VertexCollection VertexCollection;
typedef edm::View<reco::Candidate> CandidateView;

void beginStream(edm::StreamID) override {}
void produce(edm::Event &, const edm::EventSetup &) override;
Expand All @@ -42,14 +43,14 @@ class DeepBoostedJetTagInfoProducer : public edm::stream::EDProducer<> {
void fillParticleFeatures(DeepBoostedJetFeatures &fts, const reco::Jet &jet);
void fillSVFeatures(DeepBoostedJetFeatures &fts, const reco::Jet &jet);

const bool has_puppi_weighted_daughters_;
const double jet_radius_;
const double min_jet_pt_;
const double min_pt_for_track_properties_;

edm::EDGetTokenT<edm::View<reco::Jet>> jet_token_;
edm::EDGetTokenT<VertexCollection> vtx_token_;
edm::EDGetTokenT<SVCollection> sv_token_;
edm::EDGetTokenT<CandidateView> pfcand_token_;

bool use_puppi_value_map_;
bool use_pvasq_value_map_;
Expand All @@ -60,6 +61,7 @@ class DeepBoostedJetTagInfoProducer : public edm::stream::EDProducer<> {

edm::Handle<VertexCollection> vtxs_;
edm::Handle<SVCollection> svs_;
edm::Handle<CandidateView> pfcands_;
edm::ESHandle<TransientTrackBuilder> track_builder_;
edm::Handle<edm::ValueMap<float>> puppi_value_map_;
edm::Handle<edm::ValueMap<int>> pvasq_value_map_;
Expand Down Expand Up @@ -103,13 +105,13 @@ const std::vector<std::string> DeepBoostedJetTagInfoProducer::sv_features_{
};

DeepBoostedJetTagInfoProducer::DeepBoostedJetTagInfoProducer(const edm::ParameterSet &iConfig)
: has_puppi_weighted_daughters_(iConfig.getParameter<bool>("has_puppi_weighted_daughters")),
jet_radius_(iConfig.getParameter<double>("jet_radius")),
: jet_radius_(iConfig.getParameter<double>("jet_radius")),
min_jet_pt_(iConfig.getParameter<double>("min_jet_pt")),
min_pt_for_track_properties_(iConfig.getParameter<double>("min_pt_for_track_properties")),
jet_token_(consumes<edm::View<reco::Jet>>(iConfig.getParameter<edm::InputTag>("jets"))),
vtx_token_(consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
sv_token_(consumes<SVCollection>(iConfig.getParameter<edm::InputTag>("secondary_vertices"))),
pfcand_token_(consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("pf_candidates"))),
use_puppi_value_map_(false),
use_pvasq_value_map_(false) {
const auto &puppi_value_map_tag = iConfig.getParameter<edm::InputTag>("puppi_value_map");
Expand All @@ -133,12 +135,12 @@ DeepBoostedJetTagInfoProducer::~DeepBoostedJetTagInfoProducer() {}
void DeepBoostedJetTagInfoProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
// pfDeepBoostedJetTagInfos
edm::ParameterSetDescription desc;
desc.add<bool>("has_puppi_weighted_daughters", true);
desc.add<double>("jet_radius", 0.8);
desc.add<double>("min_jet_pt", 150);
desc.add<double>("min_pt_for_track_properties", -1);
desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
desc.add<edm::InputTag>("secondary_vertices", edm::InputTag("inclusiveCandidateSecondaryVertices"));
desc.add<edm::InputTag>("pf_candidates", edm::InputTag("particleFlow"));
desc.add<edm::InputTag>("jets", edm::InputTag("ak8PFJetsPuppi"));
desc.add<edm::InputTag>("puppi_value_map", edm::InputTag("puppi"));
desc.add<edm::InputTag>("vertex_associator", edm::InputTag("primaryVertexAssociation", "original"));
Expand All @@ -162,6 +164,8 @@ void DeepBoostedJetTagInfoProducer::produce(edm::Event &iEvent, const edm::Event

iEvent.getByToken(sv_token_, svs_);

iEvent.getByToken(pfcand_token_, pfcands_);

iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", track_builder_);

if (use_puppi_value_map_) {
Expand Down Expand Up @@ -252,20 +256,14 @@ void DeepBoostedJetTagInfoProducer::fillParticleFeatures(DeepBoostedJetFeatures
// remove particles w/ extremely low puppi weights
if ((puppiWgt(cand)) < 0.01)
continue;
daughters.push_back(cand);
}
// sort by (Puppi-weighted) pt
if (!has_puppi_weighted_daughters_) {
std::sort(daughters.begin(),
daughters.end(),
[&puppi_wgt_cache](const reco::CandidatePtr &a, const reco::CandidatePtr &b) {
return puppi_wgt_cache.at(a.key()) * a->pt() > puppi_wgt_cache.at(b.key()) * b->pt();
});
} else {
std::sort(daughters.begin(), daughters.end(), [](const reco::CandidatePtr &a, const reco::CandidatePtr &b) {
return a->pt() > b->pt();
});
// get the original reco/packed candidate not scaled by the puppi weight
daughters.push_back(pfcands_->ptrAt(cand.key()));
}
// sort by Puppi-weighted pt
std::sort(
daughters.begin(), daughters.end(), [&puppi_wgt_cache](const reco::CandidatePtr &a, const reco::CandidatePtr &b) {
return puppi_wgt_cache.at(a.key()) * a->pt() > puppi_wgt_cache.at(b.key()) * b->pt();
});

// reserve space
for (const auto &name : particle_features_) {
Expand All @@ -281,12 +279,8 @@ void DeepBoostedJetTagInfoProducer::fillParticleFeatures(DeepBoostedJetFeatures
const auto *packed_cand = dynamic_cast<const pat::PackedCandidate *>(&(*cand));
const auto *reco_cand = dynamic_cast<const reco::PFCandidate *>(&(*cand));

auto puppiP4 = cand->p4();
auto puppiP4 = puppi_wgt_cache.at(cand.key()) * cand->p4();
if (packed_cand) {
if (!has_puppi_weighted_daughters_) {
puppiP4 *= puppi_wgt_cache.at(cand.key());
}

float hcal_fraction = 0.;
if (packed_cand->pdgId() == 1 || packed_cand->pdgId() == 130) {
hcal_fraction = packed_cand->hcalFraction();
Expand Down