CMS 3D CMS Logo

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

#include <PFJetFilter.h>

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

Public Member Functions

 PFJetFilter (const edm::ParameterSet &)
 
 ~PFJetFilter () 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

void beginJob () override
 
void endJob () override
 
bool filter (edm::Event &, const edm::EventSetup &) override
 
double resolution (double, bool)
 
double response (double, bool)
 

Private Attributes

PFBenchmarkAlgoalgo_
 
double deltaEt_max
 
double deltaEt_min
 
double deltaR_max
 
double deltaR_min
 
unsigned int entry
 
double eta_max
 
double eta_min
 
double genPt_cut
 
edm::EDGetTokenT< edm::View< reco::Candidate > > inputRecoLabel_
 
edm::EDGetTokenT< edm::View< reco::Candidate > > inputTruthLabel_
 
double recPt_cut
 
bool verbose
 

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 13 of file PFJetFilter.h.

Constructor & Destructor Documentation

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

Definition at line 16 of file PFJetFilter.cc.

References mps_splice::entry, MuonTCMETValueMapProducer_cff::eta_max, MuonME0RecHits_cfi::eta_min, and edm::ParameterSet::getParameter().

16  {
17  inputTruthLabel_ = consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("InputTruthLabel"));
18  inputRecoLabel_ = consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("InputRecoLabel"));
19 
20  recPt_cut = iConfig.getParameter<double>("recPt");
21  genPt_cut = iConfig.getParameter<double>("genPt");
22 
23  eta_min = iConfig.getParameter<double>("minEta");
24  eta_max = iConfig.getParameter<double>("maxEta");
25 
26  deltaR_min = iConfig.getParameter<double>("deltaRMin");
27  deltaR_max = iConfig.getParameter<double>("deltaRMax");
28 
29  deltaEt_min = iConfig.getParameter<double>("minDeltaEt");
30  deltaEt_max = iConfig.getParameter<double>("maxDeltaEt");
31 
32  verbose = iConfig.getParameter<bool>("verbose");
33 
34  entry = 0;
35 }
double deltaEt_min
Definition: PFJetFilter.h:28
T getParameter(std::string const &) const
double recPt_cut
Definition: PFJetFilter.h:26
double deltaR_max
Definition: PFJetFilter.h:31
edm::EDGetTokenT< edm::View< reco::Candidate > > inputTruthLabel_
Definition: PFJetFilter.h:34
double deltaR_min
Definition: PFJetFilter.h:30
double genPt_cut
Definition: PFJetFilter.h:27
double eta_min
Definition: PFJetFilter.h:32
double deltaEt_max
Definition: PFJetFilter.h:29
edm::EDGetTokenT< edm::View< reco::Candidate > > inputRecoLabel_
Definition: PFJetFilter.h:35
unsigned int entry
Definition: PFJetFilter.h:37
double eta_max
Definition: PFJetFilter.h:33
PFJetFilter::~PFJetFilter ( )
override

Definition at line 37 of file PFJetFilter.cc.

37 {}

Member Function Documentation

void PFJetFilter::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDFilter.

Definition at line 39 of file PFJetFilter.cc.

39 {}
void PFJetFilter::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDFilter.

Definition at line 41 of file PFJetFilter.cc.

41 {}
bool PFJetFilter::filter ( edm::Event iEvent,
const edm::EventSetup iESetup 
)
overrideprivate

Definition at line 43 of file PFJetFilter.cc.

References funct::abs(), gather_cfg::cout, reco::deltaR(), promptTrackCountingComputer_cfi::deltaRmin, mps_splice::entry, reco::Candidate::eta(), MuonTCMETValueMapProducer_cff::eta_max, MuonME0RecHits_cfi::eta_min, edm::Event::getByToken(), mps_fire::i, trackingPlots::other, reco::Candidate::phi(), edm::Handle< T >::product(), EnergyCorrector::pt, reco::Candidate::pt(), ptmin, and fftjetproducer_cfi::resolution.

43  {
44  // Typedefs to use views
45  typedef edm::View<reco::Candidate> candidateCollection;
46  typedef edm::View<reco::Candidate> candidateCollection;
47 
48  const candidateCollection *truth_candidates;
49  const candidateCollection *reco_candidates;
50 
51  // ==========================================================
52  // Retrieve!
53  // ==========================================================
54 
55  // Get Truth Candidates (GenCandidates, GenJets, etc.)
57  bool isGen = iEvent.getByToken(inputTruthLabel_, truth_hnd);
58  if (!isGen) {
59  std::cout << "Warning : no Gen jets in input !" << std::endl;
60  return false;
61  }
62 
63  truth_candidates = truth_hnd.product();
64 
65  // Get Reco Candidates (PFlow, CaloJet, etc.)
67  bool isReco = iEvent.getByToken(inputRecoLabel_, reco_hnd);
68  if (!isReco) {
69  std::cout << "Warning : no Reco jets in input !" << std::endl;
70  return false;
71  }
72  reco_candidates = reco_hnd.product();
73  if (!truth_candidates || !reco_candidates)
74  return false;
75 
76  bool pass = false;
77 
78  for (unsigned int i = 0; i < reco_candidates->size(); i++) {
79  const reco::Candidate *particle = &(*reco_candidates)[i];
80 
81  // This protection is meant at not being used !
82  assert(particle != nullptr);
83 
84  double rec_pt = particle->pt();
85  double rec_eta = particle->eta();
86  double rec_phi = particle->phi();
87 
88  // skip PFjets with pt < recPt_cut GeV
89  if (rec_pt < recPt_cut)
90  continue;
91 
92  // skip PFjets with eta > maxEta_cut or eta < minEta_cut
93  if (fabs(rec_eta) > eta_max)
94  continue;
95  if (fabs(rec_eta) < eta_min)
96  continue;
97 
98  bool Barrel = false;
99  bool Endcap = false;
100  if (std::abs(rec_eta) < 1.4)
101  Barrel = true;
102  if (std::abs(rec_eta) > 1.4 && std::abs(rec_eta) < 2.6)
103  Endcap = true;
104 
105  // Keep only jets in the barrel or the endcaps, within the tracker
106  // acceptance
107  if (!Barrel && !Endcap)
108  continue;
109 
110  // Find the closets recJet
111  double deltaRmin = 999.;
112  double ptmin = 0.;
113  for (unsigned int j = 0; j < reco_candidates->size(); j++) {
114  if (i == j)
115  continue;
116  const reco::Candidate *other = &(*reco_candidates)[j];
117  double deltaR = algo_->deltaR(particle, other);
118  if (deltaR < deltaRmin && other->pt() > 0.25 * particle->pt() && other->pt() > recPt_cut) {
119  deltaRmin = deltaR;
120  ptmin = other->pt();
121  }
122  if (deltaRmin < deltaR_min)
123  break;
124  }
125  if (deltaRmin < deltaR_min)
126  continue;
127 
128  // Find the closest genJet.
129  const reco::Candidate *gen_particle = algo_->matchByDeltaR(particle, truth_candidates);
130 
131  // Check there is a genJet associated to the recoJet
132  if (gen_particle == nullptr)
133  continue;
134 
135  // check deltaR is small enough
136  double deltaR = algo_->deltaR(particle, gen_particle);
137  if (deltaR > deltaR_max)
138  continue;
139 
140  // double true_E = gen_particle->p();
141  double true_pt = gen_particle->pt();
142  double true_eta = gen_particle->eta();
143  double true_phi = gen_particle->phi();
144 
145  // skip PFjets with pt < genPt_cut GeV
146  if (true_pt < genPt_cut)
147  continue;
148 
149  // Find the pT residual
150  double resPt = (rec_pt - true_pt) / true_pt;
151  double sigma = resolution(true_pt, Barrel);
152  double avera = response(true_pt, Barrel);
153  double nSig = (resPt - avera) / sigma;
154 
155  if (nSig > deltaEt_max || nSig < deltaEt_min) {
156  /* */
157  if (verbose)
158  std::cout << "Entry " << entry << " resPt = " << resPt << " sigma/avera/nSig = " << sigma << "/" << avera << "/"
159  << nSig << " pT (T/R) " << true_pt << "/" << rec_pt << " Eta (T/R) " << true_eta << "/" << rec_eta
160  << " Phi (T/R) " << true_phi << "/" << rec_phi << " deltaRMin/ptmin " << deltaRmin << "/" << ptmin
161  << std::endl;
162  /* */
163  pass = true;
164  }
165 
166  if (pass)
167  break;
168  }
169 
170  entry++;
171  return pass;
172 }
double deltaEt_min
Definition: PFJetFilter.h:28
PFBenchmarkAlgo * algo_
Definition: PFJetFilter.h:40
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
double recPt_cut
Definition: PFJetFilter.h:26
double deltaR_max
Definition: PFJetFilter.h:31
double response(double, bool)
Definition: PFJetFilter.cc:183
edm::EDGetTokenT< edm::View< reco::Candidate > > inputTruthLabel_
Definition: PFJetFilter.h:34
deltaRmin
maximum deltaR of track to jet.
static double deltaR(const T *, const U *)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR_min
Definition: PFJetFilter.h:30
static const Collection::value_type * matchByDeltaR(const T *, const Collection *)
virtual double eta() const =0
momentum pseudorapidity
T const * product() const
Definition: Handle.h:74
virtual double pt() const =0
transverse momentum
double genPt_cut
Definition: PFJetFilter.h:27
double eta_min
Definition: PFJetFilter.h:32
double ptmin
Definition: HydjetWrapper.h:90
double deltaEt_max
Definition: PFJetFilter.h:29
edm::EDGetTokenT< edm::View< reco::Candidate > > inputRecoLabel_
Definition: PFJetFilter.h:35
unsigned int entry
Definition: PFJetFilter.h:37
virtual double phi() const =0
momentum azimuthal angle
double resolution(double, bool)
Definition: PFJetFilter.cc:174
double eta_max
Definition: PFJetFilter.h:33
double PFJetFilter::resolution ( double  pt,
bool  barrel 
)
private

Definition at line 174 of file PFJetFilter.cc.

References p1, p2, EnergyCorrector::pt, and mathSSE::sqrt().

174  {
175  double p0 = barrel ? 1.19200e-02 : 8.45341e-03;
176  double p1 = barrel ? 1.06138e+00 : 7.96855e-01;
177  double p2 = barrel ? -2.05929e+00 : -3.12076e-01;
178 
179  double resp = p0 + p1 / sqrt(pt) + p2 / pt;
180  return resp;
181 }
T sqrt(T t)
Definition: SSEVec.h:18
double p2[4]
Definition: TauolaWrapper.h:90
double p1[4]
Definition: TauolaWrapper.h:89
double PFJetFilter::response ( double  pt,
bool  barrel 
)
private

Definition at line 183 of file PFJetFilter.cc.

References JetChargeProducer_cfi::exp, p1, and p2.

183  {
184  double p0 = barrel ? 1.09906E-1 : 6.91939E+1;
185  double p1 = barrel ? -1.61443E-1 : -6.92733E+1;
186  double p2 = barrel ? 3.45489E+3 : 1.58207E+6;
187 
188  double resp = p0 + p1 * std::exp(-pt / p2);
189  return resp;
190 }
double p2[4]
Definition: TauolaWrapper.h:90
double p1[4]
Definition: TauolaWrapper.h:89

Member Data Documentation

PFBenchmarkAlgo* PFJetFilter::algo_
private

Definition at line 40 of file PFJetFilter.h.

double PFJetFilter::deltaEt_max
private

Definition at line 29 of file PFJetFilter.h.

double PFJetFilter::deltaEt_min
private

Definition at line 28 of file PFJetFilter.h.

double PFJetFilter::deltaR_max
private

Definition at line 31 of file PFJetFilter.h.

double PFJetFilter::deltaR_min
private

Definition at line 30 of file PFJetFilter.h.

unsigned int PFJetFilter::entry
private

Definition at line 37 of file PFJetFilter.h.

double PFJetFilter::eta_max
private

Definition at line 33 of file PFJetFilter.h.

double PFJetFilter::eta_min
private

Definition at line 32 of file PFJetFilter.h.

double PFJetFilter::genPt_cut
private

Definition at line 27 of file PFJetFilter.h.

edm::EDGetTokenT<edm::View<reco::Candidate> > PFJetFilter::inputRecoLabel_
private

Definition at line 35 of file PFJetFilter.h.

edm::EDGetTokenT<edm::View<reco::Candidate> > PFJetFilter::inputTruthLabel_
private

Definition at line 34 of file PFJetFilter.h.

double PFJetFilter::recPt_cut
private

Definition at line 26 of file PFJetFilter.h.

bool PFJetFilter::verbose
private

Definition at line 38 of file PFJetFilter.h.