CMS 3D CMS Logo

V0EventSelector.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <vector>
12 
14 public:
15  explicit V0EventSelector(const edm::ParameterSet&);
16  ~V0EventSelector() override = default;
17 
18  bool filter(edm::Event&, const edm::EventSetup&) override;
19  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
20 
21 private:
23  const unsigned int minNumCandidates_;
24  const double massMin_;
25  const double massMax_;
26 };
27 
29  : vccToken_{consumes<reco::VertexCompositeCandidateCollection>(
30  iConfig.getParameter<edm::InputTag>("vertexCompositeCandidates"))},
31  minNumCandidates_{iConfig.getParameter<unsigned int>("minCandidates")},
32  massMin_{iConfig.getParameter<double>("massMin")},
33  massMax_{iConfig.getParameter<double>("massMax")} {
34  produces<reco::VertexCompositeCandidateCollection>();
35 }
36 
39  iEvent.getByToken(vccToken_, vccHandle);
40  auto filteredVCC = std::make_unique<reco::VertexCompositeCandidateCollection>();
41 
42  // early return if the input collection is empty
43  if (!vccHandle.isValid()) {
44  iEvent.put(std::move(filteredVCC));
45  return false;
46  }
47 
48  for (const auto& vcc : *vccHandle) {
49  if (vcc.mass() >= massMin_ && vcc.mass() <= massMax_) {
50  filteredVCC->push_back(vcc);
51  }
52  }
53 
54  bool pass = filteredVCC->size() >= minNumCandidates_;
55  iEvent.put(std::move(filteredVCC));
56 
57  return pass;
58 }
59 
62  desc.add<edm::InputTag>("vertexCompositeCandidates", edm::InputTag("generalV0Candidates:Kshort"));
63  desc.add<unsigned int>("minCandidates", 1)->setComment("Minimum number of candidates required");
64  desc.add<double>("massMin", 0.0)->setComment("Minimum mass in GeV");
65  desc.add<double>("massMax", std::numeric_limits<double>::max())->setComment("Maximum mass in GeV");
66  descriptions.addWithDefaultLabel(desc);
67 }
68 
69 // Define this module as a plug-in
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
~V0EventSelector() override=default
int iEvent
Definition: GenABIO.cc:224
bool filter(edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
V0EventSelector(const edm::ParameterSet &)
const edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > vccToken_
const double massMin_
const unsigned int minNumCandidates_
bool isValid() const
Definition: HandleBase.h:70
const double massMax_
def move(src, dest)
Definition: eostools.py:511