CMS 3D CMS Logo

PtMinPFCandidateSelectorDefinition.h
Go to the documentation of this file.
1 #ifndef CommonTools_ParticleFlow_PtMinPFCandidateSelectorDefinition
2 #define CommonTools_ParticleFlow_PtMinPFCandidateSelectorDefinition
3 
8 
9 namespace pf2pat {
10 
12 
13  public:
15  ptMin_( cfg.getParameter< double >( "ptMin" ) ) { }
16 
17 
18  void select( const HandleToCollection & hc,
19  const edm::EventBase & e,
20  const edm::EventSetup& s
21  ) {
22  selected_.clear();
23 
24  assert( hc.isValid() );
25 
26 
27  unsigned key=0;
28  for( collection::const_iterator pfc = hc->begin();
29  pfc != hc->end(); ++pfc, ++key) {
30 
31  if( pfc->pt() > ptMin_ ) {
32  selected_.push_back( reco::PFCandidate(*pfc) );
33  reco::PFCandidatePtr ptrToMother( hc, key );
34  selected_.back().setSourceCandidatePtr( ptrToMother );
35 
36  }
37  }
38  }
39 
40 /* const container& selected() const {return selected_;} */
41 
42  private:
43  double ptMin_;
44  };
45 
46 }
47 
48 #endif
void select(const HandleToCollection &hc, const edm::EventBase &e, const edm::EventSetup &s)
bool isValid() const
Definition: HandleBase.h:74
PtMinPFCandidateSelectorDefinition(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
susybsm::HSCParticleCollection hc
Definition: classes.h:25