CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
VBFGenJetFilter Class Reference

#include <VBFGenJetFilter.h>

Inheritance diagram for VBFGenJetFilter:
edm::EDFilter edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

bool filter (edm::Event &, const edm::EventSetup &) override
 
 VBFGenJetFilter (const edm::ParameterSet &)
 
 ~VBFGenJetFilter () override
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDFilter () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

int charge (const int &Id)
 
std::vector< const reco::GenJet * > filterGenJets (const std::vector< reco::GenJet > *jets)
 
std::vector< const reco::GenParticle * > filterGenLeptons (const std::vector< reco::GenParticle > *particles)
 
std::vector< HepMC::GenParticle * > getNu (const HepMC::GenEvent *particles)
 
std::vector< HepMC::GenParticle * > getSt3 (const HepMC::GenEvent *particles)
 
std::vector< HepMC::GenParticle * > getVisibleDecayProducts (HepMC::GenParticle *particle)
 
double nuMET (std::vector< HepMC::GenParticle * > vNu)
 
void printGenVector (std::vector< HepMC::GenParticle * > vec)
 

Private Attributes

double deltaRJetLep
 
double etaMax
 
double etaMin
 
bool leadJetsNoLepMass
 
edm::EDGetTokenT< reco::GenJetCollectionm_inputTag_GenJetCollection
 
edm::EDGetTokenT< reco::GenParticleCollectionm_inputTag_GenParticleCollection
 
double maxDeltaEta
 
double maxDeltaPhi
 
double maxInvMass
 
double maxLeadingJetsInvMass
 
double minDeltaEta
 
double minDeltaPhi
 
double minInvMass
 
double minLeadingJetsInvMass
 
bool oppositeHemisphere
 
double ptMin
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 29 of file VBFGenJetFilter.h.

Constructor & Destructor Documentation

VBFGenJetFilter::VBFGenJetFilter ( const edm::ParameterSet iConfig)
explicit

Definition at line 18 of file VBFGenJetFilter.cc.

References edm::ParameterSet::getUntrackedParameter(), leadJetsNoLepMass, m_inputTag_GenJetCollection, and m_inputTag_GenParticleCollection.

18  :
19 oppositeHemisphere(iConfig.getUntrackedParameter<bool> ("oppositeHemisphere",false)),
20 leadJetsNoLepMass (iConfig.getUntrackedParameter<bool> ("leadJetsNoLepMass", false)),
21 ptMin (iConfig.getUntrackedParameter<double>("minPt", 20)),
22 etaMin (iConfig.getUntrackedParameter<double>("minEta", -5.0)),
23 etaMax (iConfig.getUntrackedParameter<double>("maxEta", 5.0)),
24 minInvMass (iConfig.getUntrackedParameter<double>("minInvMass", 0.0)),
25 maxInvMass (iConfig.getUntrackedParameter<double>("maxInvMass", 99999.0)),
26 minDeltaPhi (iConfig.getUntrackedParameter<double>("minDeltaPhi", -1.0)),
27 maxDeltaPhi (iConfig.getUntrackedParameter<double>("maxDeltaPhi", 99999.0)),
28 minDeltaEta (iConfig.getUntrackedParameter<double>("minDeltaEta", -1.0)),
29 maxDeltaEta (iConfig.getUntrackedParameter<double>("maxDeltaEta", 99999.0)),
30 minLeadingJetsInvMass (iConfig.getUntrackedParameter<double>("minLeadingJetsInvMass", 0.0)),
31 maxLeadingJetsInvMass (iConfig.getUntrackedParameter<double>("maxLeadingJetsInvMass", 99999.0)),
32 deltaRJetLep (iConfig.getUntrackedParameter<double>("deltaRJetLep", 0.3))
33 {
34 
35  m_inputTag_GenJetCollection = consumes<reco::GenJetCollection>(iConfig.getUntrackedParameter<edm::InputTag>("inputTag_GenJetCollection",edm::InputTag("ak5GenJetsNoNu")));
36  if (leadJetsNoLepMass) m_inputTag_GenParticleCollection = consumes<reco::GenParticleCollection>(iConfig.getUntrackedParameter<edm::InputTag>("genParticles",edm::InputTag("genParticles")));
37 
38 }
T getUntrackedParameter(std::string const &, T const &) const
double minLeadingJetsInvMass
edm::EDGetTokenT< reco::GenParticleCollection > m_inputTag_GenParticleCollection
double maxLeadingJetsInvMass
edm::EDGetTokenT< reco::GenJetCollection > m_inputTag_GenJetCollection
VBFGenJetFilter::~VBFGenJetFilter ( )
override

Definition at line 40 of file VBFGenJetFilter.cc.

40  {
41 
42 }

Member Function Documentation

int VBFGenJetFilter::charge ( const int &  Id)
private
bool VBFGenJetFilter::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 84 of file VBFGenJetFilter.cc.

References a, b, DEFINE_FWK_MODULE, reco::deltaPhi(), reco::deltaR2(), deltaRJetLep, particleFlow_cfi::dEta, particleFlow_cfi::dPhi, filterGenJets(), filterGenLeptons(), ttbarCategorization_cff::genJets, GenHFHadronMatcher_cfi::genParticles, edm::Event::getByToken(), mps_fire::i, boostedTaus_cff::jetIdx, leadJetsNoLepMass, m_inputTag_GenJetCollection, m_inputTag_GenParticleCollection, maxDeltaEta, maxDeltaPhi, maxInvMass, maxLeadingJetsInvMass, minLeadingJetsInvMass, oppositeHemisphere, p4, reco::LeafCandidate::p4(), and edm::Handle< T >::product().

85 {
86  using namespace edm;
87 
88  Handle< vector<reco::GenJet> > handleGenJets;
89  iEvent.getByToken(m_inputTag_GenJetCollection, handleGenJets);
90  const vector<reco::GenJet>* genJets = handleGenJets.product();
91 
92  // Getting filtered generator jets
93  vector<const reco::GenJet*> filGenJets = filterGenJets(genJets);
94 
95  // If we do not find at least 2 jets veto the event
96  if(filGenJets.size()<2){return false;}
97 
98 
99  // Testing dijet mass
100  if(leadJetsNoLepMass) {
101 
102  Handle<reco::GenParticleCollection> genParticelesCollection;
103  iEvent.getByToken(m_inputTag_GenParticleCollection, genParticelesCollection);
104  const vector<reco::GenParticle>* genParticles = genParticelesCollection.product();
105 
106  // Getting filtered generator muons
107  vector<const reco::GenParticle*> filGenLep = filterGenLeptons(genParticles);
108 
109  // Getting p4 of jet with no lepton
110  vector<math::XYZTLorentzVector> genJetsWithoutLeptonsP4;
111  unsigned int jetIdx = 0;
112 
113  while(genJetsWithoutLeptonsP4.size()<2 && jetIdx < filGenJets.size()) {
114  bool jetWhitoutLep = true;
115 
116  const math::XYZTLorentzVector & p4J= (filGenJets[jetIdx])->p4();
117  for(unsigned int i = 0; i < filGenLep.size() && jetWhitoutLep; ++i) {
118  if(reco::deltaR2((filGenLep[i])->p4(), p4J) < deltaRJetLep*deltaRJetLep)
119  jetWhitoutLep = false;
120  }
121  if (jetWhitoutLep) genJetsWithoutLeptonsP4.push_back(p4J);
122  ++jetIdx;
123  }
124 
125  // Checking the invariant mass of the leading jets
126  if (genJetsWithoutLeptonsP4.size() < 2) return false;
127  float invMassLeadingJet = (genJetsWithoutLeptonsP4[0] + genJetsWithoutLeptonsP4[1]).M();
128  if ( invMassLeadingJet > minLeadingJetsInvMass && invMassLeadingJet < maxLeadingJetsInvMass) return true;
129  else return false;
130  }
131 
132 
133  for(unsigned a=0; a<filGenJets.size(); a++){
134  for(unsigned b=a+1; b<filGenJets.size(); b++){
135 
136  const reco::GenJet* pA = filGenJets[a];
137  const reco::GenJet* pB = filGenJets[b];
138 
139  // Getting the dijet vector
140  math::XYZTLorentzVector diJet = pA->p4() + pB->p4();
141 
142  // Testing opposite hemispheres
143  double dijetProd = pA->p4().eta()*pB->p4().eta();
144  if(oppositeHemisphere && dijetProd>=0){continue;}
145 
146  // Testing dijet mass
147  double invMass = diJet.mass();
148  if(invMass<=minInvMass || invMass>maxInvMass){continue;}
149 
150  // Testing dijet delta eta
151  double dEta = fabs(pA->p4().eta()-pB->p4().eta());
152  if(dEta<=minDeltaEta || dEta>maxDeltaEta){continue;}
153 
154  // Testing dijet delta phi
155  double dPhi = fabs(reco::deltaPhi(pA->p4().phi(),pB->p4().phi()));
156  if(dPhi<=minDeltaPhi || dPhi>maxDeltaPhi){continue;}
157 
158  return true;
159  }
160  }
161 
162  return false;
163 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
double minLeadingJetsInvMass
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< const reco::GenJet * > filterGenJets(const std::vector< reco::GenJet > *jets)
edm::EDGetTokenT< reco::GenParticleCollection > m_inputTag_GenParticleCollection
std::vector< const reco::GenParticle * > filterGenLeptons(const std::vector< reco::GenParticle > *particles)
double maxLeadingJetsInvMass
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double p4[4]
Definition: TauolaWrapper.h:92
Jets made from MC generator particles.
Definition: GenJet.h:25
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
T const * product() const
Definition: Handle.h:74
double b
Definition: hdecay.h:120
HLT enums.
double a
Definition: hdecay.h:121
edm::EDGetTokenT< reco::GenJetCollection > m_inputTag_GenJetCollection
vector< const reco::GenJet * > VBFGenJetFilter::filterGenJets ( const std::vector< reco::GenJet > *  jets)
private

Definition at line 65 of file VBFGenJetFilter.cc.

References etaMax, etaMin, mps_fire::i, MillePedeFileConverter_cfg::out, reco::LeafCandidate::p4(), and ptMin.

Referenced by filter().

65  {
66 
67  vector<const reco::GenJet*> out;
68 
69  for(unsigned i=0; i<jets->size(); i++){
70 
71  const reco::GenJet* j = &((*jets)[i]);
72 
73  if(j->p4().pt() >ptMin && j->p4().eta()>etaMin && j->p4().eta()<etaMax)
74  {
75  out.push_back(j);
76  }
77  }
78 
79  return out;
80 }
vector< PseudoJet > jets
Jets made from MC generator particles.
Definition: GenJet.h:25
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
vector< const reco::GenParticle * > VBFGenJetFilter::filterGenLeptons ( const std::vector< reco::GenParticle > *  particles)
private

Definition at line 45 of file VBFGenJetFilter.cc.

References funct::abs(), MillePedeFileConverter_cfg::out, and AlCaHLTBitMon_ParallelJobs::p.

Referenced by filter().

45  {
46  vector<const reco::GenParticle*> out;
47 
48 
49 
50  for(const auto & p : *particles){
51 
52  int absPdgId = std::abs(p.pdgId());
53 
54  if(((absPdgId == 11) || (absPdgId == 13) || (absPdgId == 15)) && p.isHardProcess()) {
55  out.push_back(&p);
56  }
57 
58 
59  }
60  return out;
61 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector<HepMC::GenParticle*> VBFGenJetFilter::getNu ( const HepMC::GenEvent particles)
private
std::vector<HepMC::GenParticle*> VBFGenJetFilter::getSt3 ( const HepMC::GenEvent particles)
private
std::vector<HepMC::GenParticle*> VBFGenJetFilter::getVisibleDecayProducts ( HepMC::GenParticle *  particle)
private
double VBFGenJetFilter::nuMET ( std::vector< HepMC::GenParticle * >  vNu)
private
void VBFGenJetFilter::printGenVector ( std::vector< HepMC::GenParticle * >  vec)
private

Member Data Documentation

double VBFGenJetFilter::deltaRJetLep
private

Definition at line 68 of file VBFGenJetFilter.h.

Referenced by filter().

double VBFGenJetFilter::etaMax
private

Definition at line 59 of file VBFGenJetFilter.h.

Referenced by filterGenJets().

double VBFGenJetFilter::etaMin
private

Definition at line 58 of file VBFGenJetFilter.h.

Referenced by filterGenJets().

bool VBFGenJetFilter::leadJetsNoLepMass
private

Definition at line 56 of file VBFGenJetFilter.h.

Referenced by filter(), and VBFGenJetFilter().

edm::EDGetTokenT< reco::GenJetCollection > VBFGenJetFilter::m_inputTag_GenJetCollection
private

Definition at line 71 of file VBFGenJetFilter.h.

Referenced by filter(), and VBFGenJetFilter().

edm::EDGetTokenT< reco::GenParticleCollection > VBFGenJetFilter::m_inputTag_GenParticleCollection
private

Definition at line 72 of file VBFGenJetFilter.h.

Referenced by filter(), and VBFGenJetFilter().

double VBFGenJetFilter::maxDeltaEta
private

Definition at line 65 of file VBFGenJetFilter.h.

Referenced by filter().

double VBFGenJetFilter::maxDeltaPhi
private

Definition at line 63 of file VBFGenJetFilter.h.

Referenced by filter().

double VBFGenJetFilter::maxInvMass
private

Definition at line 61 of file VBFGenJetFilter.h.

Referenced by filter().

double VBFGenJetFilter::maxLeadingJetsInvMass
private

Definition at line 67 of file VBFGenJetFilter.h.

Referenced by filter().

double VBFGenJetFilter::minDeltaEta
private

Definition at line 64 of file VBFGenJetFilter.h.

double VBFGenJetFilter::minDeltaPhi
private

Definition at line 62 of file VBFGenJetFilter.h.

double VBFGenJetFilter::minInvMass
private

Definition at line 60 of file VBFGenJetFilter.h.

double VBFGenJetFilter::minLeadingJetsInvMass
private

Definition at line 66 of file VBFGenJetFilter.h.

Referenced by filter().

bool VBFGenJetFilter::oppositeHemisphere
private

Definition at line 55 of file VBFGenJetFilter.h.

Referenced by filter().

double VBFGenJetFilter::ptMin
private

Definition at line 57 of file VBFGenJetFilter.h.

Referenced by filterGenJets().