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
Updated vbf filter v2 #20394
Updated vbf filter v2 #20394
Changes from 3 commits
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -16,6 +16,7 @@ using namespace std; | |
|
||
VBFGenJetFilter::VBFGenJetFilter(const edm::ParameterSet& iConfig) : | ||
oppositeHemisphere(iConfig.getUntrackedParameter<bool> ("oppositeHemisphere",false)), | ||
leadJetsNoLepMass (iConfig.getUntrackedParameter<bool> ("leadJetsNoLepMass", false)), | ||
ptMin (iConfig.getUntrackedParameter<double>("minPt", 20)), | ||
etaMin (iConfig.getUntrackedParameter<double>("minEta", -5.0)), | ||
etaMax (iConfig.getUntrackedParameter<double>("maxEta", 5.0)), | ||
|
@@ -24,17 +25,41 @@ maxInvMass (iConfig.getUntrackedParameter<double>("maxInvMass", 9999 | |
minDeltaPhi (iConfig.getUntrackedParameter<double>("minDeltaPhi", -1.0)), | ||
maxDeltaPhi (iConfig.getUntrackedParameter<double>("maxDeltaPhi", 99999.0)), | ||
minDeltaEta (iConfig.getUntrackedParameter<double>("minDeltaEta", -1.0)), | ||
maxDeltaEta (iConfig.getUntrackedParameter<double>("maxDeltaEta", 99999.0)) | ||
maxDeltaEta (iConfig.getUntrackedParameter<double>("maxDeltaEta", 99999.0)), | ||
minLeadingJetsInvMass (iConfig.getUntrackedParameter<double>("minLeadingJetsInvMass", 0.0)), | ||
maxLeadingJetsInvMass (iConfig.getUntrackedParameter<double>("maxLeadingJetsInvMass", 99999.0)), | ||
deltaRJetLep (iConfig.getUntrackedParameter<double>("deltaRJetLep", 0.3)) | ||
{ | ||
|
||
m_inputTag_GenJetCollection = consumes<reco::GenJetCollection>(iConfig.getUntrackedParameter<edm::InputTag>("inputTag_GenJetCollection",edm::InputTag("ak5GenJetsNoNu"))); | ||
m_inputTag_GenParticleCollection = consumes<reco::GenParticleCollection>(iConfig.getUntrackedParameter<edm::InputTag>("genParticles",edm::InputTag("genParticles"))); | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. this is only consumed if leadJetsNoLepMass is true - please add that conditional here. |
||
|
||
} | ||
|
||
VBFGenJetFilter::~VBFGenJetFilter(){ | ||
|
||
} | ||
|
||
|
||
vector<const reco::GenParticle*> VBFGenJetFilter::filterGenLeptons(const vector<reco::GenParticle>* particles){ | ||
vector<const reco::GenParticle*> out; | ||
|
||
for(unsigned i=0; i<particles->size(); i++){ | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. hi @gmandorl - please use a range based loop here (for const & p: *particles) { |
||
|
||
const reco::GenParticle* p = &((*particles)[i]); | ||
int absPdgId = abs(p->pdgId()); | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. please use std::abs here |
||
|
||
if(((absPdgId == 11) || (absPdgId == 13) || (absPdgId == 15)) && p->isHardProcess()) { | ||
out.push_back(p); | ||
} | ||
|
||
|
||
} | ||
return out; | ||
} | ||
|
||
|
||
|
||
vector<const reco::GenJet*> VBFGenJetFilter::filterGenJets(const vector<reco::GenJet>* jets){ | ||
|
||
vector<const reco::GenJet*> out; | ||
|
@@ -68,6 +93,42 @@ bool VBFGenJetFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) | |
// If we do not find at least 2 jets veto the event | ||
if(filGenJets.size()<2){return false;} | ||
|
||
|
||
// Testing dijet mass | ||
if(leadJetsNoLepMass) { | ||
|
||
Handle<reco::GenParticleCollection> genParticelesCollection; | ||
iEvent.getByToken(m_inputTag_GenParticleCollection, genParticelesCollection); | ||
const vector<reco::GenParticle>* genParticles = genParticelesCollection.product(); | ||
|
||
// Getting filtered generator muons | ||
vector<const reco::GenParticle*> filGenLep = filterGenLeptons(genParticles); | ||
|
||
// Getting p4 of jet with no lepton | ||
vector<math::XYZTLorentzVector> genJetsWithoutLeptonsP4; | ||
unsigned int jetIdx = 0; | ||
|
||
while(genJetsWithoutLeptonsP4.size()<2 && jetIdx < filGenJets.size()) { | ||
bool jetWhitoutLep = true; | ||
math::XYZTLorentzVector p4J= (filGenJets[jetIdx])->p4(); | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. could this be a const &? |
||
for(unsigned int i = 0; i < filGenLep.size() && jetWhitoutLep; ++i) { | ||
float dPhi = reco::deltaPhi((filGenLep[i])->p4().phi(), p4J.phi()); | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. is this code not just reco::deltaR2? |
||
float dEta = (filGenLep[i])->p4().eta()-p4J.eta(); | ||
if(dPhi*dPhi + dEta*dEta < deltaRJetLep*deltaRJetLep) | ||
jetWhitoutLep = false; | ||
} | ||
if (jetWhitoutLep) genJetsWithoutLeptonsP4.push_back(p4J); | ||
++jetIdx; | ||
} | ||
|
||
// Checking the invariant mass of the leading jets | ||
if (genJetsWithoutLeptonsP4.size() < 2) return false; | ||
float invMassLeadingJet = (genJetsWithoutLeptonsP4[0] + genJetsWithoutLeptonsP4[1]).M(); | ||
if ( invMassLeadingJet > minLeadingJetsInvMass && invMassLeadingJet < maxLeadingJetsInvMass) return true; | ||
else return false; | ||
} | ||
|
||
|
||
for(unsigned a=0; a<filGenJets.size(); a++){ | ||
for(unsigned b=a+1; b<filGenJets.size(); b++){ | ||
|
||
|
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
this should probably be False
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
not true, should be True...
Does it work if leadJetsNoLepMass and the other parameters are not specified?
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
If leadJetsNoLepMass is true but no other parameters are specified the code works.
In that case, the events pass the selection whenever there are at least two jets with DR > 0.3 from all the charged hard leptons.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
If leadJetsNoLepMass is not specified, it is set as false