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

Add isolated charged hadron flag to RECO+AOD #18324

Merged
merged 3 commits into from Apr 20, 2017
Merged
Show file tree
Hide file tree
Changes from 2 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
Expand Up @@ -44,7 +44,8 @@
'keep *_particleFlow_photons_*',
'keep *_trackerDrivenElectronSeeds_preid_*',
'keep *_particleFlowPtrs_*_*',
'keep *_particleFlowTmpPtrs_*_*'
'keep *_particleFlowTmpPtrs_*_*',
'keep *_chargedHadronPFTrackIsolation_*_*'
)
)
# RECO content
Expand Down Expand Up @@ -88,7 +89,8 @@
'keep *_particleFlow_muons_*',
'keep *_trackerDrivenElectronSeeds_preid_*',
'keep *_particleFlowPtrs_*_*',
'keep *_particleFlowTmpPtrs_*_*'
'keep *_particleFlowTmpPtrs_*_*',
'keep *_chargedHadronPFTrackIsolation_*_*'
)
)

Expand Down Expand Up @@ -130,7 +132,8 @@
'keep recoPhotonCores_pfPhotonTranslator_*_*',
'keep recoConversions_pfPhotonTranslator_*_*',
'keep *_particleFlowPtrs_*_*',
'keep *_particleFlowTmpPtrs_*_*'
'keep *_particleFlowTmpPtrs_*_*',
'keep *_chargedHadronPFTrackIsolation_*_*'
)
)

Expand Down
Expand Up @@ -19,6 +19,7 @@
from CommonTools.ParticleFlow.pfParticleSelection_cff import *

from RecoEgamma.EgammaIsolationAlgos.particleBasedIsoProducer_cff import *
from RecoParticleFlow.PFProducer.chargedHadronPFTrackIsolation_cfi import *

from RecoJets.JetProducers.fixedGridRhoProducerFastjet_cfi import *
fixedGridRhoFastjetAllTmp = fixedGridRhoFastjetAll.clone(pfCandidatesTag = cms.InputTag("particleFlowTmp"))
Expand All @@ -36,7 +37,7 @@
particleFlowEGammaFinal*
pfParticleSelectionSequence )

particleFlowLinks = cms.Sequence( particleFlow*particleFlowPtrs*particleBasedIsolationSequence)
particleFlowLinks = cms.Sequence( particleFlow*particleFlowPtrs*chargedHadronPFTrackIsolation*particleBasedIsolationSequence)

from RecoParticleFlow.PFTracking.hgcalTrackCollection_cfi import *
from RecoParticleFlow.PFProducer.simPFProducer_cfi import *
Expand Down
@@ -0,0 +1,101 @@
/*
* ChargedHadronPFTrackIsolationProducer
*
* Author: Andreas Hinzmann
*
* Associates PF isolation flag to charged hadron candidates
*
*/

#include "FWCore/Framework/interface/stream/EDProducer.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"

#include "FWCore/Utilities/interface/Exception.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"

#include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
#include "DataFormats/Common/interface/ValueMap.h"
#include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
#include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"

namespace edm {
class ConfigurationDescriptions;
}
Copy link
Contributor

Choose a reason for hiding this comment

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

this forward declaration can be removed (there is an include for this class earlier in this file)


class ChargedHadronPFTrackIsolationProducer : public edm::stream::EDProducer<>
Copy link
Contributor

Choose a reason for hiding this comment

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

this module can be made edm::global

Copy link
Contributor Author

Choose a reason for hiding this comment

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

made edm::global and changed the produce method to "void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& es) const" as defined in global/EDProducerBase.h.

{
public:

explicit ChargedHadronPFTrackIsolationProducer(const edm::ParameterSet& cfg);
~ChargedHadronPFTrackIsolationProducer() {}
void produce(edm::Event& evt, const edm::EventSetup& es);
Copy link
Contributor

Choose a reason for hiding this comment

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

please add override at the end

static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);

private:
// input collection
edm::InputTag srcCandidates_;
edm::EDGetTokenT<edm::View<reco::PFCandidate> > Candidates_token;
Copy link
Contributor

@slava77 slava77 Apr 18, 2017

Choose a reason for hiding this comment

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

please keep the same style of data members: start with lower case, end with an underscore

double minTrackPt_;
double minRawCaloEnergy_;

};

ChargedHadronPFTrackIsolationProducer::ChargedHadronPFTrackIsolationProducer(const edm::ParameterSet& cfg)
{
srcCandidates_ = cfg.getParameter<edm::InputTag>("src");
Candidates_token = consumes<edm::View<reco::PFCandidate> >(srcCandidates_);
minTrackPt_ = cfg.getParameter<double>("minTrackPt");
minRawCaloEnergy_ = cfg.getParameter<double>("minRawCaloEnergy");

produces<edm::ValueMap<bool> >();
}

void ChargedHadronPFTrackIsolationProducer::produce(edm::Event& evt, const edm::EventSetup& es)
{
// get a view of our Candidates via the base candidates
typedef edm::View<reco::PFCandidate> PFCandidateView;
edm::Handle<PFCandidateView> Candidates;
evt.getByToken(Candidates_token, Candidates);
Copy link
Contributor

Choose a reason for hiding this comment

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

start local variables and data members with lower case letter. Capitalization is reserved for class names.


std::vector<bool> values;

for( PFCandidateView::const_iterator c = Candidates->begin(); c!=Candidates->end(); ++c) {
Copy link
Contributor

Choose a reason for hiding this comment

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

for range loop seems more appropriate for (auto cons& c : *candidates) {

// Check that there is only one track in the block.
unsigned int nTracks = 0;
if ((c->particleId()==1) && (c->pt()>minTrackPt_) && ((c->rawEcalEnergy()+c->rawHcalEnergy())>minRawCaloEnergy_)) {
const reco::PFCandidate::ElementsInBlocks& theElements = c->elementsInBlocks();
if( theElements.empty() ) continue;
const reco::PFBlockRef blockRef = theElements[0].first;
const edm::OwnVector<reco::PFBlockElement>& elements = blockRef->elements();
for(unsigned int iEle=0; iEle<elements.size(); iEle++) { // Find the tracks in the block
Copy link
Contributor

Choose a reason for hiding this comment

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

similarly

 for (auto const& ele : elements){
...

reco::PFBlockElement::Type type = elements[iEle].type();
if( type== reco::PFBlockElement::TRACK)
nTracks++;
}
}
values.push_back((nTracks==1));
}

std::unique_ptr<edm::ValueMap<bool> > out(new edm::ValueMap<bool>());
edm::ValueMap<bool>::Filler filler(*out);
filler.insert(Candidates,values.begin(),values.end());
filler.fill();
evt.put(std::move(out));
}

void ChargedHadronPFTrackIsolationProducer::fillDescriptions(edm::ConfigurationDescriptions & descriptions) {
edm::ParameterSetDescription desc;
desc.add<edm::InputTag>("src", edm::InputTag("particleFlow"));
desc.add<double>("minTrackPt", 1);
desc.add<double>("minRawCaloEnergy", 0.5);
descriptions.add("chargedHadronPFTrackIsolation", desc);
}

#include "FWCore/Framework/interface/MakerMacros.h"

DEFINE_FWK_MODULE(ChargedHadronPFTrackIsolationProducer);
@@ -0,0 +1,9 @@
import FWCore.ParameterSet.Config as cms


chargedHadronPFTrackIsolation = cms.EDProducer(
Copy link
Contributor

Choose a reason for hiding this comment

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

this file is redundant and can be removed.
It is already provided by the framework with the call descriptions.add("chargedHadronPFTrackIsolation", desc); available in ChargedHadronPFTrackIsolationProducer

"ChargedHadronPFTrackIsolationProducer",
src = cms.InputTag("particleFlow"),
minTrackPt = cms.double(1),
minRawCaloEnergy = cms.double(0.5),
)