CMS 3D CMS Logo

NearbyCandCountComputer.cc
Go to the documentation of this file.
1 //
2 //
3 
19 
23 
27 
29 
31  public:
32  explicit NearbyCandCountComputer(const edm::ParameterSet & iConfig);
33  virtual ~NearbyCandCountComputer() ;
34 
35  virtual void produce(edm::Event & iEvent, const edm::EventSetup & iSetup) override;
36 
37  private:
40  double deltaR2_;
41  StringCutObjectSelector<reco::Candidate,true> objCut_; // lazy parsing, to allow cutting on variables not in reco::Candidate class
43 };
44 
46  probesToken_(consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("probes"))),
47  objectsToken_(consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("objects"))),
48  deltaR2_(std::pow(iConfig.getParameter<double>("deltaR"), 2)),
49  objCut_(iConfig.existsAs<std::string>("objectSelection") ? iConfig.getParameter<std::string>("objectSelection") : "", true),
50  pairCut_(iConfig.existsAs<std::string>("pairSelection") ? iConfig.getParameter<std::string>("pairSelection") : "", true)
51 {
52  produces<edm::ValueMap<float> >();
53 }
54 
55 
57 {
58 }
59 
60 void
62  using namespace edm;
63 
64  // read input
66  iEvent.getByToken(probesToken_, probes);
67  iEvent.getByToken(objectsToken_, objects);
68 
69  // prepare vector for output
70  std::vector<float> values;
71 
72  // fill
73  View<reco::Candidate>::const_iterator probe, endprobes = probes->end();
74  View<reco::Candidate>::const_iterator object, beginobjects = objects->begin(), endobjects = objects->end();
75  for (probe = probes->begin(); probe != endprobes; ++probe) {
76  float count = 0;
77  for (object = beginobjects; object != endobjects; ++object) {
78  if ((deltaR2(*probe, *object) < deltaR2_) &&
79  objCut_(*object) &&
80  pairCut_(pat::DiObjectProxy(*probe, *object))) {
81  count++;
82  }
83  }
84  values.push_back(count);
85  }
86 
87  // convert into ValueMap and store
88  auto valMap = std::make_unique<ValueMap<float>>();
90  filler.insert(probes, values.begin(), values.end());
91  filler.fill();
92  iEvent.put(std::move(valMap));
93 }
94 
95 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
StringCutObjectSelector< pat::DiObjectProxy, true > pairCut_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::View< reco::Candidate > > objectsToken_
StringCutObjectSelector< reco::Candidate, true > objCut_
edm::EDGetTokenT< edm::View< reco::Candidate > > probesToken_
Count candidates near to another candidate, write result in ValueMap.
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
NearbyCandCountComputer(const edm::ParameterSet &iConfig)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:510