Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:53

0001 /**
0002   \class    pat::PATSingleVertexSelector PATSingleVertexSelector.h "PhysicsTools/PatAlgos/plugins/PATSingleVertexSelector.h"
0003   \brief    Produces a list containing a single vertex selected by some criteria
0004 
0005 
0006   \author   Giovanni Petrucciani
0007   \version  $Id: PATSingleVertexSelector.h,v 1.5 2011/06/15 11:47:25 friis Exp $
0008 */
0009 
0010 #include "DataFormats/Common/interface/View.h"
0011 #include "DataFormats/Candidate/interface/VertexCompositeCandidate.h"
0012 #include <DataFormats/BeamSpot/interface/BeamSpot.h>
0013 #include "FWCore/Utilities/interface/transform.h"
0014 #include "FWCore/Framework/interface/stream/EDFilter.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/Utilities/interface/InputTag.h"
0018 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0019 #include "DataFormats/VertexReco/interface/Vertex.h"
0020 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0021 #include "DataFormats/Candidate/interface/Candidate.h"
0022 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0023 
0024 #include <algorithm>
0025 
0026 namespace pat {
0027 
0028   class PATSingleVertexSelector : public edm::stream::EDFilter<> {
0029   public:
0030     explicit PATSingleVertexSelector(const edm::ParameterSet &iConfig);
0031     ~PATSingleVertexSelector() override;
0032 
0033     bool filter(edm::Event &iEvent, const edm::EventSetup &iSetup) override;
0034 
0035   private:
0036     enum Mode { First, NearestToCand, FromCand, FromBeamSpot };
0037     typedef StringCutObjectSelector<reco::Vertex> VtxSel;
0038     typedef StringCutObjectSelector<reco::Candidate> CandSel;
0039 
0040     Mode parseMode(const std::string &name) const;
0041 
0042     std::unique_ptr<std::vector<reco::Vertex>> filter_(Mode mode,
0043                                                        const edm::Event &iEvent,
0044                                                        const edm::EventSetup &iSetup);
0045     bool hasMode_(Mode mode) const;
0046     // configurables
0047     std::vector<Mode> modes_;  // mode + optional fallbacks
0048     edm::EDGetTokenT<std::vector<reco::Vertex>> verticesToken_;
0049     std::vector<edm::EDGetTokenT<edm::View<reco::Candidate>>> candidatesToken_;
0050     const VtxSel vtxPreselection_;
0051     const CandSel candPreselection_;
0052     edm::EDGetTokenT<reco::BeamSpot> beamSpotToken_;
0053     // transient data. meaningful while 'filter()' is on the stack
0054     std::vector<reco::VertexRef> selVtxs_;
0055     reco::CandidatePtr bestCand_;
0056 
0057     // flag to enable/disable EDFilter functionality:
0058     // if set to false, PATSingleVertexSelector selects the "one" event vertex,
0059     // but does not reject any events
0060     bool doFilterEvents_;
0061   };
0062 
0063 }  // namespace pat
0064 
0065 using pat::PATSingleVertexSelector;
0066 
0067 PATSingleVertexSelector::Mode PATSingleVertexSelector::parseMode(const std::string &mode) const {
0068   if (mode == "firstVertex") {
0069     return First;
0070   } else if (mode == "nearestToCandidate") {
0071     return NearestToCand;
0072   } else if (mode == "fromCandidate") {
0073     return FromCand;
0074   } else if (mode == "beamSpot") {
0075     return FromBeamSpot;
0076   } else {
0077     throw cms::Exception("Configuration")
0078         << "PATSingleVertexSelector: Mode '" << mode << "' not recognized or not supported.\n";
0079   }
0080 }
0081 
0082 PATSingleVertexSelector::PATSingleVertexSelector(const edm::ParameterSet &iConfig)
0083     : vtxPreselection_(iConfig.existsAs<std::string>("vertexPreselection")
0084                            ? iConfig.getParameter<std::string>("vertexPreselection")
0085                            : std::string(" 1 == 1 ")),
0086       candPreselection_(iConfig.existsAs<std::string>("candidatePreselection")
0087                             ? iConfig.getParameter<std::string>("candidatePreselection")
0088                             : std::string(" 1 == 1 ")),
0089       doFilterEvents_(false) {
0090   using namespace std;
0091 
0092   modes_.push_back(parseMode(iConfig.getParameter<std::string>("mode")));
0093   if (iConfig.exists("fallbacks")) {
0094     vector<string> modes = iConfig.getParameter<vector<string>>("fallbacks");
0095     for (vector<string>::const_iterator it = modes.begin(), ed = modes.end(); it != ed; ++it) {
0096       modes_.push_back(parseMode(*it));
0097     }
0098   }
0099   if (hasMode_(First) || hasMode_(NearestToCand)) {
0100     verticesToken_ = consumes<vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("vertices"));
0101   }
0102   if (hasMode_(NearestToCand) || hasMode_(FromCand)) {
0103     candidatesToken_ =
0104         edm::vector_transform(iConfig.getParameter<vector<edm::InputTag>>("candidates"),
0105                               [this](edm::InputTag const &tag) { return consumes<edm::View<reco::Candidate>>(tag); });
0106   }
0107   if (hasMode_(FromBeamSpot)) {
0108     beamSpotToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"));
0109   }
0110 
0111   if (iConfig.exists("filter"))
0112     doFilterEvents_ = iConfig.getParameter<bool>("filter");
0113 
0114   produces<vector<reco::Vertex>>();
0115 }
0116 
0117 PATSingleVertexSelector::~PATSingleVertexSelector() {}
0118 
0119 bool PATSingleVertexSelector::hasMode_(Mode mode) const {
0120   return (std::find(modes_.begin(), modes_.end(), mode) != modes_.end());
0121 }
0122 
0123 bool PATSingleVertexSelector::filter(edm::Event &iEvent, const edm::EventSetup &iSetup) {
0124   using namespace edm;
0125   using namespace std;
0126 
0127   // Clear
0128   selVtxs_.clear();
0129   bestCand_ = reco::CandidatePtr();
0130 
0131   // Gather data from the Event
0132   // -- vertex data --
0133   if (hasMode_(First) || hasMode_(NearestToCand)) {
0134     Handle<vector<reco::Vertex>> vertices;
0135     iEvent.getByToken(verticesToken_, vertices);
0136     for (vector<reco::Vertex>::const_iterator itv = vertices->begin(), edv = vertices->end(); itv != edv; ++itv) {
0137       if (!(vtxPreselection_(*itv)))
0138         continue;
0139       selVtxs_.push_back(reco::VertexRef(vertices, std::distance(vertices->begin(), itv)));
0140     }
0141   }
0142   // -- candidate data --
0143   if (hasMode_(NearestToCand) || hasMode_(FromCand)) {
0144     vector<pair<double, reco::CandidatePtr>> cands;
0145     for (vector<edm::EDGetTokenT<edm::View<reco::Candidate>>>::const_iterator itt = candidatesToken_.begin(),
0146                                                                               edt = candidatesToken_.end();
0147          itt != edt;
0148          ++itt) {
0149       Handle<View<reco::Candidate>> theseCands;
0150       iEvent.getByToken(*itt, theseCands);
0151       for (View<reco::Candidate>::const_iterator itc = theseCands->begin(), edc = theseCands->end(); itc != edc;
0152            ++itc) {
0153         if (!(candPreselection_(*itc)))
0154           continue;
0155         cands.push_back(pair<double, reco::CandidatePtr>(
0156             -itc->pt(), reco::CandidatePtr(theseCands, std::distance(theseCands->begin(), itc))));
0157       }
0158     }
0159     if (!cands.empty())
0160       bestCand_ = cands.front().second;
0161   }
0162 
0163   bool passes = false;
0164   std::unique_ptr<vector<reco::Vertex>> result;
0165   // Run main mode + possible fallback modes
0166   for (std::vector<Mode>::const_iterator itm = modes_.begin(), endm = modes_.end(); itm != endm; ++itm) {
0167     result = filter_(*itm, iEvent, iSetup);
0168     // Check if we got any vertices.  If so, take them.
0169     if (!result->empty()) {
0170       passes = true;
0171       break;
0172     }
0173   }
0174   iEvent.put(std::move(result));
0175   // Check if we want to apply the EDFilter
0176   if (doFilterEvents_)
0177     return passes;
0178   else
0179     return true;
0180 }
0181 
0182 std::unique_ptr<std::vector<reco::Vertex>> PATSingleVertexSelector::filter_(Mode mode,
0183                                                                             const edm::Event &iEvent,
0184                                                                             const edm::EventSetup &iSetup) {
0185   using namespace edm;
0186   using namespace std;
0187   auto result = std::make_unique<std::vector<reco::Vertex>>();
0188   switch (mode) {
0189     case First: {
0190       if (selVtxs_.empty())
0191         return result;
0192       result->push_back(*selVtxs_.front());
0193       return result;
0194     }
0195     case FromCand: {
0196       if (bestCand_.isNull())
0197         return result;
0198       reco::Vertex vtx;
0199       auto const &bestCandDeref = *bestCand_;
0200       if (typeid(bestCandDeref) == typeid(reco::VertexCompositeCandidate)) {
0201         vtx = reco::Vertex(bestCand_->vertex(),
0202                            bestCand_->vertexCovariance(),
0203                            bestCand_->vertexChi2(),
0204                            bestCand_->vertexNdof(),
0205                            bestCand_->numberOfDaughters());
0206       } else {
0207         vtx = reco::Vertex(bestCand_->vertex(), reco::Vertex::Error(), 0, 0, 0);
0208       }
0209       result->push_back(vtx);
0210       return result;
0211     }
0212     case NearestToCand: {
0213       if (selVtxs_.empty() || (bestCand_.isNull()))
0214         return result;
0215       reco::VertexRef which;
0216       float dzmin = 9999.0;
0217       for (auto itv = selVtxs_.begin(), edv = selVtxs_.end(); itv != edv; ++itv) {
0218         float dz = std::abs((*itv)->z() - bestCand_->vz());
0219         if (dz < dzmin) {
0220           dzmin = dz;
0221           which = *itv;
0222         }
0223       }
0224       if (which.isNonnull())  // actually it should not happen, but better safe than sorry
0225         result->push_back(*which);
0226       return result;
0227     }
0228     case FromBeamSpot: {
0229       Handle<reco::BeamSpot> beamSpot;
0230       iEvent.getByToken(beamSpotToken_, beamSpot);
0231       reco::Vertex bs(beamSpot->position(), beamSpot->covariance3D(), 0, 0, 0);
0232       result->push_back(bs);
0233       return result;
0234     }
0235     default:
0236       // Return an empty vector signifying no vertices found.
0237       return result;
0238   }
0239 }
0240 
0241 #include "FWCore/Framework/interface/MakerMacros.h"
0242 DEFINE_FWK_MODULE(PATSingleVertexSelector);