Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:09

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 // make a selector for this selection
0014 class PVSelector : public Selector<edm::EventBase> {
0015 public:
0016   PVSelector() {}
0017 
0018   PVSelector(edm::ParameterSet const& params) : pvSrc_(params.getParameter<edm::InputTag>("pvSrc")), pvSel_(params) {
0019     push_back("NPV", params.getParameter<int>("NPV"));
0020     set("NPV");
0021     retInternal_ = getBitTemplate();
0022     indexNPV_ = index_type(&bits_, "NPV");
0023   }
0024 
0025 #ifndef __GCCXML__
0026   PVSelector(edm::ParameterSet const& params, edm::ConsumesCollector&& iC) : PVSelector(params) {
0027     pvSrcToken_ = iC.consumes<std::vector<reco::Vertex> >(pvSrc_);
0028   }
0029 #endif
0030 
0031   bool operator()(edm::EventBase const& event, pat::strbitset& ret) override {
0032     ret.set(false);
0033     event.getByLabel(pvSrc_, h_primVtx);
0034 
0035     // check if there is a good primary vertex
0036 
0037     if (h_primVtx->empty())
0038       return false;
0039 
0040     // Loop over PV's and count those that pass
0041     int npv = 0;
0042     int _ntotal = 0;
0043     mvSelPvs.clear();
0044     for (std::vector<reco::Vertex>::const_iterator ibegin = h_primVtx->begin(), iend = h_primVtx->end(), i = ibegin;
0045          i != iend;
0046          ++i) {
0047       reco::Vertex const& pv = *i;
0048       bool ipass = pvSel_(pv);
0049       if (ipass) {
0050         ++npv;
0051         mvSelPvs.push_back(edm::Ptr<reco::Vertex>(h_primVtx, _ntotal));
0052       }
0053       ++_ntotal;
0054     }
0055 
0056     // cache npv
0057     mNpv = npv;
0058 
0059     // Set the strbitset
0060     if (npv >= cut(indexNPV_, int()) || ignoreCut(indexNPV_)) {
0061       passCut(ret, indexNPV_);
0062     }
0063 
0064     // Check if there is anything to ignore
0065     setIgnored(ret);
0066 
0067     // Return status
0068     bool pass = (bool)ret;
0069     return pass;
0070   }
0071 
0072   using EventSelector::operator();
0073 
0074   edm::Handle<std::vector<reco::Vertex> > const& vertices() const { return h_primVtx; }
0075 
0076   // get NPV from the last check
0077   int GetNpv(void) { return mNpv; }
0078 
0079   std::vector<edm::Ptr<reco::Vertex> > const& GetSelectedPvs() const { return mvSelPvs; }
0080 
0081 private:
0082   edm::InputTag pvSrc_;
0083 #ifndef __GCCXML__
0084   edm::EDGetTokenT<std::vector<reco::Vertex> > pvSrcToken_;
0085 #endif
0086   PVObjectSelector pvSel_;
0087   edm::Handle<std::vector<reco::Vertex> > h_primVtx;
0088   std::vector<edm::Ptr<reco::Vertex> > mvSelPvs;  // selected vertices
0089   index_type indexNPV_;
0090   int mNpv;  // cache number of PVs
0091 };
0092 
0093 #endif