CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
IsolatedPFCandidateSelectorDefinition.h
Go to the documentation of this file.
1 #ifndef CommonTools_ParticleFlow_IsolatedPFCandidateSelectorDefinition
2 #define CommonTools_ParticleFlow_IsolatedPFCandidateSelectorDefinition
3 
11 
12 namespace pf2pat {
13 
15 
16  public:
18 
20  isolationValueMapChargedTokens_(edm::vector_transform(cfg.getParameter< std::vector<edm::InputTag> >("isolationValueMapsCharged"), [&](edm::InputTag const & tag){return iC.consumes<IsoMap>(tag);}) ),
21  isolationValueMapNeutralTokens_(edm::vector_transform(cfg.getParameter< std::vector<edm::InputTag> >("isolationValueMapsNeutral"), [&](edm::InputTag const & tag){return iC.consumes<IsoMap>(tag);}) ),
22  doDeltaBetaCorrection_(cfg.getParameter<bool>("doDeltaBetaCorrection")),
23  deltaBetaIsolationValueMapToken_(iC.mayConsume<IsoMap>(cfg.getParameter< edm::InputTag >("deltaBetaIsolationValueMap") ) ),
24  deltaBetaFactor_(cfg.getParameter<double>("deltaBetaFactor")),
25  isRelative_(cfg.getParameter<bool>("isRelative")),
26  isolationCut_(cfg.getParameter<double>("isolationCut")) {}
27 
28 
29 
30  void select( const HandleToCollection & hc,
31  const edm::Event & e,
32  const edm::EventSetup& s) {
33  selected_.clear();
34 
35 
36  // read all charged isolation value maps
37  std::vector< edm::Handle<IsoMap> >
38  isoMapsCharged(isolationValueMapChargedTokens_.size());
39  for(unsigned iMap = 0; iMap<isolationValueMapChargedTokens_.size(); ++iMap) {
40  e.getByToken(isolationValueMapChargedTokens_[iMap], isoMapsCharged[iMap]);
41  }
42 
43 
44  // read all neutral isolation value maps
45  std::vector< edm::Handle<IsoMap> >
46  isoMapsNeutral(isolationValueMapNeutralTokens_.size());
47  for(unsigned iMap = 0; iMap<isolationValueMapNeutralTokens_.size(); ++iMap) {
48  e.getByToken(isolationValueMapNeutralTokens_[iMap], isoMapsNeutral[iMap]);
49  }
50 
51  edm::Handle<IsoMap> dBetaH;
52  if(doDeltaBetaCorrection_) {
53  e.getByToken(deltaBetaIsolationValueMapToken_, dBetaH);
54  }
55 
56  unsigned key=0;
57  for( collection::const_iterator pfc = hc->begin();
58  pfc != hc->end(); ++pfc, ++key) {
59  reco::PFCandidateRef candidate(hc,key);
60 
61  bool passed = true;
62  double isoSumCharged=0.0;
63  double isoSumNeutral=0.0;
64 
65  for(unsigned iMap = 0; iMap<isoMapsCharged.size(); ++iMap) {
66  const IsoMap & isoMap = *(isoMapsCharged[iMap]);
67  double val = isoMap[candidate];
68  isoSumCharged+=val;
69  }
70 
71  for(unsigned iMap = 0; iMap<isoMapsNeutral.size(); ++iMap) {
72  const IsoMap & isoMap = *(isoMapsNeutral[iMap]);
73  double val = isoMap[candidate];
74  isoSumNeutral+=val;
75  }
76 
77 
78  if ( doDeltaBetaCorrection_ ) {
79  const IsoMap& isoMap = *dBetaH;
80  double dBetaVal = isoMap[candidate];
81  double dBetaCorIsoSumNeutral = isoSumNeutral + deltaBetaFactor_*dBetaVal;
82  isoSumNeutral = dBetaCorIsoSumNeutral>0 ? dBetaCorIsoSumNeutral : 0; //follow muon POG definition in 2012
83  }
84 
85  double isoSum=isoSumCharged+isoSumNeutral;
86 
87  if( isRelative_ ) {
88  isoSum /= candidate->pt();
89  }
90 
91  if ( isoSum>isolationCut_ ) {
92  passed = false;
93  }
94 
95  if(passed) {
96  // passed all cuts, selected
97  selected_.push_back( reco::PFCandidate(*pfc) );
98  reco::PFCandidatePtr ptrToMother( hc, key );
99  selected_.back().setSourceCandidatePtr( ptrToMother );
100  }
101  }
102  }
103 
104 
105  private:
106  std::vector<edm::EDGetTokenT<IsoMap> > isolationValueMapChargedTokens_;
107  std::vector<edm::EDGetTokenT<IsoMap> > isolationValueMapNeutralTokens_;
113  };
114 
115 }
116 
117 #endif
IsolatedPFCandidateSelectorDefinition(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
tuple cfg
Definition: looper.py:293
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
isolationCut_(cfg.getParameter< double >("isolationCut"))
tuple isoSum
===&gt; require is not the leading lepton and opposite to the leading lepton
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::vector< edm::EDGetTokenT< IsoMap > > isolationValueMapChargedTokens_
void select(const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s)
string key
FastSim: produces sample of signal events, overlayed with premixed minbias events.
std::vector< edm::EDGetTokenT< IsoMap > > isolationValueMapNeutralTokens_
string const
Definition: compareJSON.py:14
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
susybsm::HSCParticleCollection hc
Definition: classes.h:25