Skip to content

Commit

Permalink
Merge pull request #10751 from ahinzmann/partonFlavourMiniAOD75
Browse files Browse the repository at this point in the history
Switch parton flavour in MiniAOD from algorithmic to physics definition 75
  • Loading branch information
cmsbuild committed Sep 2, 2015
2 parents b58c545 + 0927982 commit 5c6a928
Show file tree
Hide file tree
Showing 7 changed files with 27 additions and 11 deletions.
16 changes: 13 additions & 3 deletions PhysicsTools/JetMCAlgos/plugins/HadronAndPartonSelector.cc
Expand Up @@ -118,7 +118,8 @@ HadronAndPartonSelector::HadronAndPartonSelector(const edm::ParameterSet& iConfi
//register your products
produces<reco::GenParticleRefVector>( "bHadrons" );
produces<reco::GenParticleRefVector>( "cHadrons" );
produces<reco::GenParticleRefVector>( "partons" );
produces<reco::GenParticleRefVector>( "algorithmicPartons" );
produces<reco::GenParticleRefVector>( "physicsPartons" );
produces<reco::GenParticleRefVector>( "leptons" );
}

Expand Down Expand Up @@ -206,6 +207,7 @@ HadronAndPartonSelector::produce(edm::Event& iEvent, const edm::EventSetup& iSet
std::auto_ptr<reco::GenParticleRefVector> bHadrons ( new reco::GenParticleRefVector );
std::auto_ptr<reco::GenParticleRefVector> cHadrons ( new reco::GenParticleRefVector );
std::auto_ptr<reco::GenParticleRefVector> partons ( new reco::GenParticleRefVector );
std::auto_ptr<reco::GenParticleRefVector> physicsPartons ( new reco::GenParticleRefVector );
std::auto_ptr<reco::GenParticleRefVector> leptons ( new reco::GenParticleRefVector );

// loop over particles and select b and c hadrons and leptons
Expand Down Expand Up @@ -249,12 +251,20 @@ HadronAndPartonSelector::produce(edm::Event& iEvent, const edm::EventSetup& iSet
}

// select partons
if ( partonMode_!="Undefined" )
if ( partonMode_!="Undefined" ) {
partonSelector_->run(particles,partons);
for(reco::GenParticleCollection::const_iterator it = particles->begin(); it != particles->end(); ++it)
{
if( !(it->status()==3 || (( partonMode_=="Pythia8" ) && (it->status()==23)))) continue;
if( !CandMCTagUtils::isParton( *it ) ) continue; // skip particle if not a parton
physicsPartons->push_back( reco::GenParticleRef( particles, it - particles->begin() ) );
}
}

iEvent.put( bHadrons, "bHadrons" );
iEvent.put( cHadrons, "cHadrons" );
iEvent.put( partons, "partons" );
iEvent.put( partons, "algorithmicPartons" );
iEvent.put( physicsPartons, "physicsPartons" );
iEvent.put( leptons, "leptons" );
}

Expand Down
Expand Up @@ -4,7 +4,8 @@
jets = cms.InputTag("ak4PFJets"),
bHadrons = cms.InputTag("selectedHadronsAndPartons","bHadrons"),
cHadrons = cms.InputTag("selectedHadronsAndPartons","cHadrons"),
partons = cms.InputTag("selectedHadronsAndPartons","partons"),
partons = cms.InputTag("selectedHadronsAndPartons","physicsPartons"),
leptons = cms.InputTag("selectedHadronsAndPartons","leptons"),
jetAlgorithm = cms.string("AntiKt"),
rParam = cms.double(0.4),
ghostRescaling = cms.double(1e-18),
Expand Down
Expand Up @@ -4,7 +4,8 @@
jets = cms.InputTag("ak5PFJets"),
bHadrons = cms.InputTag("selectedHadronsAndPartons","bHadrons"),
cHadrons = cms.InputTag("selectedHadronsAndPartons","cHadrons"),
partons = cms.InputTag("selectedHadronsAndPartons","partons"),
partons = cms.InputTag("selectedHadronsAndPartons","physicsPartons"),
leptons = cms.InputTag("selectedHadronsAndPartons","leptons"),
jetAlgorithm = cms.string("AntiKt"),
rParam = cms.double(0.5),
ghostRescaling = cms.double(1e-18),
Expand Down
Expand Up @@ -26,11 +26,12 @@
jets = cms.InputTag("ak4PFJetsCHS"),
bHadrons = cms.InputTag("patJetPartons","bHadrons"),
cHadrons = cms.InputTag("patJetPartons","cHadrons"),
partons = cms.InputTag("patJetPartons","partons"),
partons = cms.InputTag("patJetPartons","physicsPartons"),
leptons = cms.InputTag("patJetPartons","leptons"),
jetAlgorithm = cms.string("AntiKt"),
rParam = cms.double(0.4),
ghostRescaling = cms.double(1e-18),
hadronFlavourHasPriority = cms.bool(True)
hadronFlavourHasPriority = cms.bool(False)
)

# default PAT sequence for jet flavour identification
Expand Down
Expand Up @@ -75,7 +75,7 @@
# jet flavour idetification configurables
getJetMCFlavour = cms.bool(True),
useLegacyJetMCFlavour = cms.bool(False),
addJetFlavourInfo = cms.bool(False),
addJetFlavourInfo = cms.bool(True),
JetPartonMapSource = cms.InputTag("patJetFlavourAssociationLegacy"),
JetFlavourInfoSource = cms.InputTag("patJetFlavourAssociation"),
# efficiencies
Expand Down
3 changes: 2 additions & 1 deletion PhysicsTools/PatAlgos/python/slimming/miniAOD_tools.py
Expand Up @@ -299,6 +299,7 @@ def miniAOD_customizeMC(process):
process.photonMatch.src = cms.InputTag("reducedEgamma","reducedGedPhotons")
process.tauMatch.matched = "prunedGenParticles"
process.tauGenJets.GenParticles = "prunedGenParticles"
process.patJetPartons.particles = "prunedGenParticles"
process.patJetPartonMatch.matched = "prunedGenParticles"
process.patJetPartonMatch.mcStatus = [ 3, 23 ]
process.patJetGenJetMatch.matched = "slimmedGenJets"
Expand All @@ -308,7 +309,7 @@ def miniAOD_customizeMC(process):
process.patPhotons.embedGenMatch = False
process.patTaus.embedGenMatch = False
process.patJets.embedGenPartonMatch = False
#also jet flavour must be switched to ak4
#also jet flavour must be switched
process.patJetFlavourAssociation.rParam = 0.4

def miniAOD_customizeOutput(out):
Expand Down
6 changes: 4 additions & 2 deletions PhysicsTools/PatAlgos/python/tools/jetTools.py
Expand Up @@ -318,7 +318,8 @@ def toolCode(self, process):
_newPatJetFlavourAssociation.rParam=rParam
_newPatJetFlavourAssociation.bHadrons=cms.InputTag("patJetPartons"+postfix,"bHadrons")
_newPatJetFlavourAssociation.cHadrons=cms.InputTag("patJetPartons"+postfix,"cHadrons")
_newPatJetFlavourAssociation.partons=cms.InputTag("patJetPartons"+postfix,"partons")
_newPatJetFlavourAssociation.partons=cms.InputTag("patJetPartons"+postfix,"physicsPartons")
_newPatJetFlavourAssociation.leptons=cms.InputTag("patJetPartons"+postfix,"leptons")
else :
setattr(process, 'patJetFlavourAssociation'+_labelName+postfix,
patJetFlavourAssociation.clone(
Expand All @@ -327,7 +328,8 @@ def toolCode(self, process):
rParam=rParam,
bHadrons = cms.InputTag("patJetPartons"+postfix,"bHadrons"),
cHadrons = cms.InputTag("patJetPartons"+postfix,"cHadrons"),
partons = cms.InputTag("patJetPartons"+postfix,"partons")
partons = cms.InputTag("patJetPartons"+postfix,"physicsPartons"),
leptons = cms.InputTag("patJetPartons"+postfix,"leptons")
)
)
knownModules.append('patJetFlavourAssociation'+_labelName+postfix)
Expand Down

0 comments on commit 5c6a928

Please sign in to comment.