Skip to content

Commit

Permalink
Merge pull request #47015 from Dr15Jones/lowMemPuppi
Browse files Browse the repository at this point in the history
Decrease temporary memory held by PuppiProducer
  • Loading branch information
cmsbuild authored Dec 20, 2024
2 parents 666efe2 + c382271 commit a852d33
Show file tree
Hide file tree
Showing 5 changed files with 147 additions and 161 deletions.
4 changes: 2 additions & 2 deletions CommonTools/PileupAlgos/interface/PuppiAlgo.h
Original file line number Diff line number Diff line change
Expand Up @@ -9,7 +9,7 @@

class PuppiAlgo {
public:
PuppiAlgo(edm::ParameterSet &iConfig);
PuppiAlgo(const edm::ParameterSet &iConfig);
~PuppiAlgo();
static void fillDescriptionsPuppiAlgo(edm::ParameterSetDescription &desc);
//Computing Mean and RMS
Expand All @@ -19,7 +19,7 @@ class PuppiAlgo {
void computeMedRMS(const unsigned int &iAlgo);
//Get the Weight
double compute(std::vector<double> const &iVals, double iChi2) const;
const std::vector<float> &alphas() { return fPups; }
const std::vector<float> &alphas() const { return fPups; }
//Helpers
inline int etaBins() const { return fEtaMin.size(); }
inline double etaMin(int i) const { return fEtaMin[i]; }
Expand Down
69 changes: 33 additions & 36 deletions CommonTools/PileupAlgos/interface/PuppiContainer.h
Original file line number Diff line number Diff line change
Expand Up @@ -8,60 +8,57 @@
class PuppiContainer {
public:
PuppiContainer(const edm::ParameterSet &iConfig);
~PuppiContainer();
void initialize(const std::vector<RecoObj> &iRecoObjects);
void setPUProxy(double const iPUProxy) { fPUProxy = iPUProxy; }

std::vector<PuppiCandidate> const &pfParticles() const { return fPFParticles; }
std::vector<double> const &puppiWeights();
const std::vector<double> &puppiRawAlphas() { return fRawAlphas; }
const std::vector<double> &puppiAlphas() { return fVals; }
// const std::vector<double> puppiAlpha () {return fAlpha;}
const std::vector<double> &puppiAlphasMed() { return fAlphaMed; }
const std::vector<double> &puppiAlphasRMS() { return fAlphaRMS; }
struct Weights {
std::vector<double> weights;
std::vector<double> puppiRawAlphas;
std::vector<double> puppiAlphas;
std::vector<double> puppiAlphasMed;
std::vector<double> puppiAlphasRMS;
};

int puppiNAlgos() { return fNAlgos; }
Weights calculatePuppiWeights(const std::vector<RecoObj> &iRecoObjects, double iPUProxy);

protected:
double goodVar(PuppiCandidate const &iPart, std::vector<PuppiCandidate> const &iParts, int iOpt, const double iRCone);
int puppiNAlgos() const { return fPuppiAlgo.size(); }

private:
void initialize(const std::vector<RecoObj> &iRecoObjects,
std::vector<PuppiCandidate> &fPFParticles,
std::vector<PuppiCandidate> &fPFParticlesForVar,
std::vector<PuppiCandidate> &fPFParticlesForVarChargedPV) const;

double goodVar(PuppiCandidate const &iPart,
std::vector<PuppiCandidate> const &iParts,
int iOpt,
const double iRCone) const;
void getRMSAvg(int iOpt,
std::vector<PuppiCandidate> const &iConstits,
std::vector<PuppiCandidate> const &iParticles,
std::vector<PuppiCandidate> const &iChargeParticles);
void getRawAlphas(int iOpt,
std::vector<PuppiCandidate> const &iConstits,
std::vector<PuppiCandidate> const &iParticles,
std::vector<PuppiCandidate> const &iChargeParticles);
double getChi2FromdZ(double iDZ);
std::vector<PuppiCandidate> const &iChargeParticles,
std::vector<double> &oVals);
std::vector<double> getRawAlphas(int iOpt,
std::vector<PuppiCandidate> const &iConstits,
std::vector<PuppiCandidate> const &iParticles,
std::vector<PuppiCandidate> const &iChargeParticles) const;
double getChi2FromdZ(double iDZ) const;
int getPuppiId(float iPt, float iEta);
double var_within_R(int iId,
const std::vector<PuppiCandidate> &particles,
const PuppiCandidate &centre,
const double R);

bool fPuppiDiagnostics;
const std::vector<RecoObj> *fRecoParticles;
std::vector<PuppiCandidate> fPFParticles;
std::vector<PuppiCandidate> fPFParticlesForVar;
std::vector<PuppiCandidate> fPFParticlesForVarChargedPV;
std::vector<double> fWeights;
std::vector<double> fVals;
std::vector<double> fRawAlphas;
std::vector<double> fAlphaMed;
std::vector<double> fAlphaRMS;
const double R) const;

bool fApplyCHS;
bool fInvert;
bool fUseExp;
double fNeutralMinPt;
double fNeutralSlope;
double fPuppiWeightCut;
double fPtMaxPhotons;
double fEtaMaxPhotons;
double fPtMaxNeutrals;
double fPtMaxNeutralsStartSlope;
int fNAlgos;
double fPUProxy;
std::vector<PuppiAlgo> fPuppiAlgo;

bool fPuppiDiagnostics;
bool fApplyCHS;
bool fInvert;
bool fUseExp;
};
#endif
50 changes: 17 additions & 33 deletions CommonTools/PileupAlgos/plugins/PuppiProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -25,7 +25,6 @@
class PuppiProducer : public edm::stream::EDProducer<> {
public:
explicit PuppiProducer(const edm::ParameterSet&);
~PuppiProducer() override;

static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
typedef math::XYZTLorentzVector LorentzVector;
Expand All @@ -39,9 +38,7 @@ class PuppiProducer : public edm::stream::EDProducer<> {
typedef edm::Association<reco::VertexCollection> CandToVertex;

private:
virtual void beginJob();
void produce(edm::Event&, const edm::EventSetup&) override;
virtual void endJob();

edm::EDGetTokenT<CandidateView> tokenPFCandidates_;
edm::EDGetTokenT<VertexCollection> tokenVertices_;
Expand Down Expand Up @@ -87,12 +84,11 @@ class PuppiProducer : public edm::stream::EDProducer<> {
int fVtxNdofCut;
double fVtxZCut;
bool fUsePUProxyValue;
std::unique_ptr<PuppiContainer> fPuppiContainer;
std::vector<RecoObj> fRecoObjCollection;
PuppiContainer fPuppiContainer;
};

// ------------------------------------------------------------------------------------------
PuppiProducer::PuppiProducer(const edm::ParameterSet& iConfig) {
PuppiProducer::PuppiProducer(const edm::ParameterSet& iConfig) : fPuppiContainer(iConfig) {
fPuppiDiagnostics = iConfig.getParameter<bool>("puppiDiagnostics");
fPuppiNoLep = iConfig.getParameter<bool>("puppiNoLep");
fUseFromPVLooseTight = iConfig.getParameter<bool>("UseFromPVLooseTight");
Expand All @@ -113,7 +109,6 @@ PuppiProducer::PuppiProducer(const edm::ParameterSet& iConfig) {
fClonePackedCands = iConfig.getParameter<bool>("clonePackedCands");
fVtxNdofCut = iConfig.getParameter<int>("vtxNdofCut");
fVtxZCut = iConfig.getParameter<double>("vtxZCut");
fPuppiContainer = std::make_unique<PuppiContainer>(iConfig);

tokenPFCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
tokenVertices_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexName"));
Expand Down Expand Up @@ -150,8 +145,6 @@ PuppiProducer::PuppiProducer(const edm::ParameterSet& iConfig) {
}
}
// ------------------------------------------------------------------------------------------
PuppiProducer::~PuppiProducer() {}
// ------------------------------------------------------------------------------------------
void PuppiProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
// Get PFCandidate Collection
edm::Handle<CandidateView> hPFProduct;
Expand Down Expand Up @@ -180,11 +173,13 @@ void PuppiProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
}
}

std::vector<RecoObj> recoObjCollection;

PuppiContainer::Weights weightsInfo;
std::vector<double> lWeights;
if (!fUseExistingWeights) {
//Fill the reco objects
fRecoObjCollection.clear();
fRecoObjCollection.reserve(pfCol->size());
recoObjCollection.reserve(pfCol->size());
int iCand = 0;
for (auto const& aPF : *pfCol) {
RecoObj pReco;
Expand Down Expand Up @@ -347,15 +342,13 @@ void PuppiProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
}
}

fRecoObjCollection.push_back(pReco);
recoObjCollection.push_back(pReco);
iCand++;
}

fPuppiContainer->initialize(fRecoObjCollection);
fPuppiContainer->setPUProxy(puProxyValue);

//Compute the weights and get the particles
lWeights = fPuppiContainer->puppiWeights();
weightsInfo = fPuppiContainer.calculatePuppiWeights(recoObjCollection, puProxyValue);
lWeights = std::move(weightsInfo.weights);
} else {
//Use the existing weights
lWeights.reserve(pfCol->size());
Expand Down Expand Up @@ -457,8 +450,8 @@ void PuppiProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
p4PupFiller.insert(hPFProduct, puppiP4s.begin(), puppiP4s.end());
p4PupFiller.fill();

iEvent.emplace(ptokenPupOut_, lPupOut);
iEvent.emplace(ptokenP4PupOut_, p4PupOut);
iEvent.emplace(ptokenPupOut_, std::move(lPupOut));
iEvent.emplace(ptokenP4PupOut_, std::move(p4PupOut));
if (fUseExistingWeights || fClonePackedCands) {
edm::OrphanHandle<pat::PackedCandidateCollection> oh =
iEvent.emplace(ptokenPackedPuppiCandidates_, fPackedPuppiCandidates);
Expand All @@ -477,31 +470,22 @@ void PuppiProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
edm::ValueMap<reco::CandidatePtr>::Filler filler(pfMap_p);
filler.insert(hPFProduct, values.begin(), values.end());
filler.fill();
iEvent.emplace(ptokenValues_, pfMap_p);
iEvent.emplace(ptokenValues_, std::move(pfMap_p));

//////////////////////////////////////////////
if (fPuppiDiagnostics && !fUseExistingWeights) {
// all the different alphas per particle
// THE alpha per particle
std::vector<double> theAlphas(fPuppiContainer->puppiAlphas());
std::vector<double> theAlphasMed(fPuppiContainer->puppiAlphasMed());
std::vector<double> theAlphasRms(fPuppiContainer->puppiAlphasRMS());
std::vector<double> alphas(fPuppiContainer->puppiRawAlphas());
double nalgos(fPuppiContainer->puppiNAlgos());
double nalgos(fPuppiContainer.puppiNAlgos());

iEvent.emplace(ptokenRawAlphas_, alphas);
iEvent.emplace(ptokenRawAlphas_, std::move(weightsInfo.puppiRawAlphas));
iEvent.emplace(ptokenNalgos_, nalgos);
iEvent.emplace(ptokenAlphas_, theAlphas);
iEvent.emplace(ptokenAlphasMed_, theAlphasMed);
iEvent.emplace(ptokenAlphasRms_, theAlphasRms);
iEvent.emplace(ptokenAlphas_, std::move(weightsInfo.puppiAlphas));
iEvent.emplace(ptokenAlphasMed_, std::move(weightsInfo.puppiAlphasMed));
iEvent.emplace(ptokenAlphasRms_, std::move(weightsInfo.puppiAlphasRMS));
}
}

// ------------------------------------------------------------------------------------------
void PuppiProducer::beginJob() {}
// ------------------------------------------------------------------------------------------
void PuppiProducer::endJob() {}
// ------------------------------------------------------------------------------------------
void PuppiProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<bool>("puppiDiagnostics", false);
Expand Down
2 changes: 1 addition & 1 deletion CommonTools/PileupAlgos/src/PuppiAlgo.cc
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,7 @@
#include "Math/ProbFunc.h"
#include "TMath.h"

PuppiAlgo::PuppiAlgo(edm::ParameterSet &iConfig) {
PuppiAlgo::PuppiAlgo(edm::ParameterSet const &iConfig) {
fEtaMin = iConfig.getParameter<std::vector<double>>("etaMin");
fEtaMax = iConfig.getParameter<std::vector<double>>("etaMax");
fPtMin = iConfig.getParameter<std::vector<double>>("ptMin");
Expand Down
Loading

0 comments on commit a852d33

Please sign in to comment.