CMS 3D CMS Logo

ChargedHadronPFTrackIsolationProducer.cc
Go to the documentation of this file.
1 /*
2  * ChargedHadronPFTrackIsolationProducer
3  *
4  * Author: Andreas Hinzmann
5  *
6  * Associates PF isolation flag to charged hadron candidates
7  *
8  */
9 
16 
19 
25 
27 {
28 public:
29 
32  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
33  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
34 
35 private:
36  // input collection
39  double minTrackPt_;
41 
42 };
43 
45 {
47  candidates_token = consumes<edm::View<reco::PFCandidate> >(srccandidates_);
48  minTrackPt_ = cfg.getParameter<double>("minTrackPt");
49  minRawCaloEnergy_ = cfg.getParameter<double>("minRawCaloEnergy");
50 
51  produces<edm::ValueMap<bool> >();
52 }
53 
55 {
56  // get a view of our candidates via the base candidates
57  typedef edm::View<reco::PFCandidate> PFCandidateView;
59  evt.getByToken(candidates_token, candidates);
60 
61  std::vector<bool> values;
62 
63  for( auto const& c : *candidates ) {
64  // Check that there is only one track in the block.
65  unsigned int nTracks = 0;
66  if ((c.particleId()==1) && (c.pt()>minTrackPt_) && ((c.rawEcalEnergy()+c.rawHcalEnergy())>minRawCaloEnergy_)) {
67  const reco::PFCandidate::ElementsInBlocks& theElements = c.elementsInBlocks();
68  if( theElements.empty() ) continue;
69  const reco::PFBlockRef blockRef = theElements[0].first;
70  const edm::OwnVector<reco::PFBlockElement>& elements = blockRef->elements();
71  // Find the tracks in the block
72  for( auto const & ele : elements ) {
73  reco::PFBlockElement::Type type = ele.type();
74  if( type== reco::PFBlockElement::TRACK)
75  nTracks++;
76  }
77  }
78  values.push_back((nTracks==1));
79  }
80 
81  std::unique_ptr<edm::ValueMap<bool> > out(new edm::ValueMap<bool>());
83  filler.insert(candidates,values.begin(),values.end());
84  filler.fill();
85  evt.put(std::move(out));
86 }
87 
90  desc.add<edm::InputTag>("src", edm::InputTag("particleFlow"));
91  desc.add<double>("minTrackPt", 1);
92  desc.add<double>("minRawCaloEnergy", 0.5);
93  descriptions.add("chargedHadronPFTrackIsolation", desc);
94 }
95 
97 
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
const unsigned int nTracks(const reco::Vertex &sv)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
std::vector< ElementInBlock > ElementsInBlocks
Definition: PFCandidate.h:386
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
edm::EDGetTokenT< edm::View< reco::PFCandidate > > candidates_token
void add(std::string const &label, ParameterSetDescription const &psetDescription)
def move(src, dest)
Definition: eostools.py:510