CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ProbeMulteplicityProducer.cc
Go to the documentation of this file.
1 //
2 
15 
19 
23 
24 
26  public:
27  explicit ProbeMulteplicityProducer(const edm::ParameterSet & iConfig);
28  virtual ~ProbeMulteplicityProducer() ;
29 
30  virtual void produce(edm::Event & iEvent, const edm::EventSetup & iSetup);
31 
32  private:
34  StringCutObjectSelector<reco::Candidate,true> pairCut_; // lazy parsing, to allow cutting on variables not in reco::Candidate class
35  StringCutObjectSelector<reco::Candidate,true> probeCut_; // lazy parsing, to allow cutting on variables not in reco::Candidate class
36 };
37 
38 
40  pairs_(consumes<reco::CandidateView>(iConfig.getParameter<edm::InputTag>("pairs"))),
41  pairCut_(iConfig.existsAs<std::string>("pairSelection") ? iConfig.getParameter<std::string>("pairSelection") : "", true),
42  probeCut_(iConfig.existsAs<std::string>("probeSelection") ? iConfig.getParameter<std::string>("probeSelection") : "", true)
43 {
44  produces<edm::ValueMap<float> >();
45 }
46 
47 
48 
50 {
51 }
52 
53 
54 void
56  using namespace edm;
57 
58  // read input
60  iEvent.getByToken(pairs_, pairs);
61 
62  // fill
63  unsigned int i = 0;
64  std::vector<unsigned int> tagKeys;
65  std::vector<float> values;
66  View<reco::Candidate>::const_iterator pair, endpairs = pairs->end();
67  for (pair = pairs->begin(); pair != endpairs; ++pair, ++i) {
68  reco::CandidateBaseRef probeRef = pair->daughter(1)->masterClone();
69  unsigned int tagKey = pair->daughter(0)->masterClone().key();
70  unsigned int copies = 1;
71  if (pairCut_(*pair) && probeCut_(*probeRef)) {
72  for (unsigned int j = 0; j < i; ++j) if (tagKeys[j] == tagKey) copies++;
73  for (unsigned int j = 0; j < i; ++j) if (tagKeys[j] == tagKey) values[j] = copies;
74  } else {
76  copies = 0;
77  }
78  tagKeys.push_back(tagKey);
79  values.push_back(copies);
80  }
81 
82  // convert into ValueMap and store
83  std::auto_ptr<ValueMap<float> > valMap(new ValueMap<float>());
84  ValueMap<float>::Filler filler(*valMap);
85  filler.insert(pairs, values.begin(), values.end());
86  filler.fill();
87  iEvent.put(valMap);
88 }
89 
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::CandidateView > pairs_
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
StringCutObjectSelector< reco::Candidate, true > probeCut_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
size_t key() const
Definition: RefToBase.h:250
int j
Definition: DBlmapReader.cc:9
Matcher of number of reconstructed objects in the event to probe.
StringCutObjectSelector< reco::Candidate, true > pairCut_
const_iterator end() const
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
ProbeMulteplicityProducer(const edm::ParameterSet &iConfig)
virtual const CandidateBaseRef & masterClone() const =0