CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CandViewRefRandomSelector.cc
Go to the documentation of this file.
1 /*
2  * CandViewRefRandomSelector
3  *
4  * Author: Evan K. Friis (UC Davis)
5  *
6  * Takes a collection of objects inheriting from Candidates and returns up to
7  * N=<choose> candidates. The N output elements are selected randomly. If
8  * the collection contains N or fewer elements, the entire collection is
9  * returned.
10  *
11  */
12 
13 #include <vector>
14 
19 
22 
23 #include "TRandom3.h"
24 
26  public:
27  explicit CandViewRefRandomSelector(const edm::ParameterSet &pset);
28  bool filter(edm::Event&, const edm::EventSetup&);
29  private:
31  unsigned int choose_;
32  unsigned int seed_;
33  bool filter_;
34  TRandom3 randy_;
35 };
36 
38  const edm::ParameterSet &pset) {
39  src_ = pset.getParameter<edm::InputTag>("src");
40  choose_ = pset.getParameter<unsigned int>("choose");
41  filter_ = pset.getParameter<bool>("filter");
42  seed_ = pset.exists("seed") ? pset.getParameter<unsigned int>("seed") : 123;
43  randy_ = TRandom3(seed_);
44  produces<reco::CandidateBaseRefVector>();
45 }
46 
48  const edm::EventSetup& es) {
50  evt.getByLabel(src_, cands);
51  std::auto_ptr<reco::CandidateBaseRefVector> output(
52  new reco::CandidateBaseRefVector(cands));
53  // If we don't have enough elements to select, just copy what we have
54  const reco::CandidateBaseRefVector &inputRefs = cands->refVector();
55  if (inputRefs.size() <= choose_) {
56  *output = inputRefs;
57  } else {
58  for (size_t i = 0; i < inputRefs.size() && output->size() < choose_; ++i) {
59  // based on http://stackoverflow.com/questions/48087/
60  // select-a-random-n-elements-from-listt-in-c/48089#48089
61  double selectionProb =
62  (choose_ - output->size())*1.0/(inputRefs.size() - i);
63  // throw a number to see if we select this element
64  if (randy_.Rndm() < selectionProb)
65  output->push_back(inputRefs[i]);
66  }
67  }
68  size_t outputSize = output->size();
69  evt.put(output);
70  return ( !filter_ || outputSize );
71 }
72 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
size_type size() const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
bool filter(edm::Event &, const edm::EventSetup &)
void push_back(const RefToBase< T > &)
CandViewRefRandomSelector(const edm::ParameterSet &pset)