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

decrease event size by storing PFCandidates with half precision float… #29513

Merged
merged 3 commits into from May 5, 2020
Merged
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
13 changes: 12 additions & 1 deletion HLTrigger/JetMET/plugins/HLTScoutingPFProducer.cc
Expand Up @@ -39,6 +39,8 @@ Description: Producer for ScoutingPFJets from reco::PFJet objects, ScoutingVerte

#include "DataFormats/Math/interface/deltaR.h"

#include "DataFormats/Math/interface/libminifloat.h"

class HLTScoutingPFProducer : public edm::global::EDProducer<> {
public:
explicit HLTScoutingPFProducer(const edm::ParameterSet &);
Expand All @@ -60,6 +62,7 @@ class HLTScoutingPFProducer : public edm::global::EDProducer<> {
const double pfJetEtaCut;
const double pfCandidatePtCut;
const double pfCandidateEtaCut;
const int mantissaPrecision;

const bool doJetTags;
const bool doCandidates;
Expand All @@ -81,6 +84,7 @@ HLTScoutingPFProducer::HLTScoutingPFProducer(const edm::ParameterSet &iConfig)
pfJetEtaCut(iConfig.getParameter<double>("pfJetEtaCut")),
pfCandidatePtCut(iConfig.getParameter<double>("pfCandidatePtCut")),
pfCandidateEtaCut(iConfig.getParameter<double>("pfCandidateEtaCut")),
mantissaPrecision(iConfig.getParameter<int>("mantissaPrecision")),
doJetTags(iConfig.getParameter<bool>("doJetTags")),
doCandidates(iConfig.getParameter<bool>("doCandidates")),
doMet(iConfig.getParameter<bool>("doMet")) {
Expand Down Expand Up @@ -151,7 +155,13 @@ void HLTScoutingPFProducer::produce(edm::StreamID sid, edm::Event &iEvent, edm::
break;
++index_counter;
}
outPFCandidates->emplace_back(cand.pt(), cand.eta(), cand.phi(), cand.mass(), cand.pdgId(), vertex_index);

outPFCandidates->emplace_back(MiniFloatConverter::reduceMantissaToNbitsRounding(cand.pt(), mantissaPrecision),
MiniFloatConverter::reduceMantissaToNbitsRounding(cand.eta(), mantissaPrecision),
MiniFloatConverter::reduceMantissaToNbitsRounding(cand.phi(), mantissaPrecision),
MiniFloatConverter::reduceMantissaToNbitsRounding(cand.mass(), mantissaPrecision),
cand.pdgId(),
vertex_index);
}
}
}
Expand Down Expand Up @@ -253,6 +263,7 @@ void HLTScoutingPFProducer::fillDescriptions(edm::ConfigurationDescriptions &des
desc.add<double>("pfJetEtaCut", 3.0);
desc.add<double>("pfCandidatePtCut", 0.6);
desc.add<double>("pfCandidateEtaCut", 5.0);
desc.add<int>("mantissaPrecision", 23);
desc.add<bool>("doJetTags", true);
desc.add<bool>("doCandidates", true);
desc.add<bool>("doMet", true);
Expand Down