File indexing completed on 2025-01-12 23:42:02
0001 #ifndef Analysis_AnalysisFilters_interface_PVSelector_h
0002 #define Analysis_AnalysisFilters_interface_PVSelector_h
0003
0004 #include "FWCore/Common/interface/EventBase.h"
0005 #ifndef __GCCXML__
0006 #include "FWCore/Framework/interface/ConsumesCollector.h"
0007 #endif
0008 #include "DataFormats/Common/interface/Handle.h"
0009
0010 #include "PhysicsTools/SelectorUtils/interface/EventSelector.h"
0011 #include "PhysicsTools/SelectorUtils/interface/PVObjectSelector.h"
0012
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015
0016
0017 class PVSelector : public Selector<edm::EventBase> {
0018 public:
0019 PVSelector() {}
0020
0021 PVSelector(edm::ParameterSet const& params) : pvSrc_(params.getParameter<edm::InputTag>("pvSrc")), pvSel_(params) {
0022 push_back("NPV", params.getParameter<int>("NPV"));
0023 set("NPV");
0024 retInternal_ = getBitTemplate();
0025 indexNPV_ = index_type(&bits_, "NPV");
0026 }
0027
0028 #ifndef __GCCXML__
0029 PVSelector(edm::ParameterSet const& params, edm::ConsumesCollector&& iC) : PVSelector(params) {
0030 pvSrcToken_ = iC.consumes<std::vector<reco::Vertex> >(pvSrc_);
0031 }
0032 #endif
0033
0034 static edm::ParameterSetDescription getDescription() {
0035 edm::ParameterSetDescription desc = PVObjectSelector::getDescription();
0036 desc.add<edm::InputTag>("pvSrc", edm::InputTag(""));
0037 desc.add<int>("NPV", 1);
0038 return desc;
0039 }
0040
0041 bool operator()(edm::EventBase const& event, pat::strbitset& ret) override {
0042 ret.set(false);
0043 event.getByLabel(pvSrc_, h_primVtx);
0044
0045
0046
0047 if (h_primVtx->empty())
0048 return false;
0049
0050
0051 int npv = 0;
0052 int _ntotal = 0;
0053 mvSelPvs.clear();
0054 for (std::vector<reco::Vertex>::const_iterator ibegin = h_primVtx->begin(), iend = h_primVtx->end(), i = ibegin;
0055 i != iend;
0056 ++i) {
0057 reco::Vertex const& pv = *i;
0058 bool ipass = pvSel_(pv);
0059 if (ipass) {
0060 ++npv;
0061 mvSelPvs.push_back(edm::Ptr<reco::Vertex>(h_primVtx, _ntotal));
0062 }
0063 ++_ntotal;
0064 }
0065
0066
0067 mNpv = npv;
0068
0069
0070 if (npv >= cut(indexNPV_, int()) || ignoreCut(indexNPV_)) {
0071 passCut(ret, indexNPV_);
0072 }
0073
0074
0075 setIgnored(ret);
0076
0077
0078 bool pass = (bool)ret;
0079 return pass;
0080 }
0081
0082 using EventSelector::operator();
0083
0084 edm::Handle<std::vector<reco::Vertex> > const& vertices() const { return h_primVtx; }
0085
0086
0087 int GetNpv(void) { return mNpv; }
0088
0089 std::vector<edm::Ptr<reco::Vertex> > const& GetSelectedPvs() const { return mvSelPvs; }
0090
0091 private:
0092 edm::InputTag pvSrc_;
0093 #ifndef __GCCXML__
0094 edm::EDGetTokenT<std::vector<reco::Vertex> > pvSrcToken_;
0095 #endif
0096 PVObjectSelector pvSel_;
0097 edm::Handle<std::vector<reco::Vertex> > h_primVtx;
0098 std::vector<edm::Ptr<reco::Vertex> > mvSelPvs;
0099 index_type indexNPV_;
0100 int mNpv;
0101 };
0102
0103 #endif