CMS 3D CMS Logo

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