CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/PhysicsTools/PatAlgos/plugins/PATSingleVertexSelector.cc

Go to the documentation of this file.
00001 #include "PhysicsTools/PatAlgos/plugins/PATSingleVertexSelector.h"
00002 #include "DataFormats/Common/interface/View.h"
00003 #include "DataFormats/Candidate/interface/VertexCompositeCandidate.h"
00004 #include <DataFormats/BeamSpot/interface/BeamSpot.h>
00005 
00006 #include <algorithm>
00007 
00008 using pat::PATSingleVertexSelector;
00009 
00010 
00011 PATSingleVertexSelector::Mode
00012 PATSingleVertexSelector::parseMode(const std::string &mode) {
00013     if (mode == "firstVertex") {
00014         return First;
00015     } else if (mode == "nearestToCandidate") {
00016         return NearestToCand;
00017     } else if (mode == "fromCandidate") {
00018         return FromCand;
00019     } else if (mode == "beamSpot") {
00020         return FromBeamSpot;
00021     } else {
00022         throw cms::Exception("Configuration") << "PATSingleVertexSelector: Mode '" << mode << "' not recognized or not supported.\n";
00023     }
00024 }
00025 
00026 
00027 PATSingleVertexSelector::PATSingleVertexSelector(const edm::ParameterSet & iConfig) 
00028   : doFilterEvents_(false)
00029 {
00030    using namespace std;
00031 
00032    modes_.push_back( parseMode(iConfig.getParameter<std::string>("mode")) );
00033    if (iConfig.exists("fallbacks")) {
00034       vector<string> modes = iConfig.getParameter<vector<string> >("fallbacks");
00035       for (vector<string>::const_iterator it = modes.begin(), ed = modes.end(); it != ed; ++it) {
00036         modes_.push_back( parseMode(*it) );
00037       }
00038    }
00039    if (hasMode_(First) || hasMode_(NearestToCand)) {
00040         vertices_ = iConfig.getParameter<edm::InputTag>("vertices");
00041         if (iConfig.existsAs<string>("vertexPreselection")) {
00042             string presel = iConfig.getParameter<string>("vertexPreselection");
00043             if (!presel.empty()) vtxPreselection_ = auto_ptr<VtxSel>(new VtxSel(presel));
00044         }
00045    }
00046    if (hasMode_(NearestToCand) || hasMode_(FromCand)) {
00047         candidates_ = iConfig.getParameter<vector<edm::InputTag> >("candidates");
00048         if (iConfig.existsAs<string>("candidatePreselection")) {
00049             string presel = iConfig.getParameter<string>("candidatePreselection");
00050             if (!presel.empty()) candPreselection_ = auto_ptr<CandSel>(new CandSel(presel));
00051         }
00052    }
00053    if (hasMode_(FromBeamSpot)) {
00054         beamSpot_ = iConfig.getParameter<edm::InputTag>("beamSpot");
00055    }
00056 
00057    if ( iConfig.exists("filter") ) doFilterEvents_ = iConfig.getParameter<bool>("filter");
00058 
00059    produces<vector<reco::Vertex> >();
00060 }
00061 
00062 
00063 PATSingleVertexSelector::~PATSingleVertexSelector() 
00064 {
00065 }
00066 
00067 bool PATSingleVertexSelector::hasMode_(Mode mode) const {
00068     return (std::find(modes_.begin(), modes_.end(), mode) != modes_.end());
00069 }
00070 
00071 bool
00072 PATSingleVertexSelector::filter(edm::Event & iEvent, const edm::EventSetup & iSetup) {
00073     using namespace edm;
00074     using namespace std;
00075 
00076     // Clear
00077     selVtxs_.clear(); bestCand_ = 0;
00078 
00079     // Gather data from the Event
00080     // -- vertex data --
00081     if (hasMode_(First) || hasMode_(NearestToCand)) {
00082         Handle<vector<reco::Vertex> > vertices;
00083         iEvent.getByLabel(vertices_, vertices);
00084         for (vector<reco::Vertex>::const_iterator itv = vertices->begin(), edv = vertices->end(); itv != edv; ++itv) {
00085             if ((vtxPreselection_.get() != 0) && !((*vtxPreselection_)(*itv)) ) continue; 
00086             selVtxs_.push_back( &*itv );
00087         }
00088     }
00089     // -- candidate data --
00090     if (hasMode_(NearestToCand) || hasMode_(FromCand)) {
00091        vector<pair<double, const reco::Candidate *> > cands;
00092        for (vector<edm::InputTag>::const_iterator itt = candidates_.begin(), edt = candidates_.end(); itt != edt; ++itt) {
00093           Handle<View<reco::Candidate> > theseCands;
00094           iEvent.getByLabel(*itt, theseCands);
00095           for (View<reco::Candidate>::const_iterator itc = theseCands->begin(), edc = theseCands->end(); itc != edc; ++itc) {
00096             if ((candPreselection_.get() != 0) && !((*candPreselection_)(*itc))) continue;
00097             cands.push_back( pair<double, const reco::Candidate *>(-itc->pt(), &*itc) );
00098           }
00099        }
00100        if (!cands.empty()) bestCand_ = cands.front().second;
00101     }
00102 
00103     // Run main mode + possible fallback modes
00104     for (std::vector<Mode>::const_iterator itm = modes_.begin(), endm = modes_.end(); itm != endm; ++itm) {
00105         if (filter_(*itm, iEvent, iSetup)) return true;
00106     }
00107     if ( !doFilterEvents_ ) return true;
00108     return false;
00109 }
00110 
00111 bool
00112 PATSingleVertexSelector::filter_(Mode mode, edm::Event & iEvent, const edm::EventSetup & iSetup) {
00113     using namespace edm;
00114     using namespace std;
00115     switch(mode) {
00116         case First: {
00117             if (selVtxs_.empty()) return false;
00118             auto_ptr<vector<reco::Vertex> > result(new vector<reco::Vertex>(1, *selVtxs_.front()));
00119             iEvent.put(result);
00120             return true;
00121             }
00122         case FromCand: {
00123             if (bestCand_ == 0) return false;
00124             reco::Vertex vtx;
00125             if (typeid(*bestCand_) == typeid(reco::VertexCompositeCandidate)) {
00126                 vtx = reco::Vertex(bestCand_->vertex(), bestCand_->vertexCovariance(), 
00127                         bestCand_->vertexChi2(), bestCand_->vertexNdof(), bestCand_->numberOfDaughters() );
00128             } else {
00129                 vtx = reco::Vertex(bestCand_->vertex(), reco::Vertex::Error(), 0, 0, 0);
00130             }
00131             auto_ptr<vector<reco::Vertex> > result(new vector<reco::Vertex>(1, vtx));
00132             iEvent.put(result);
00133             return true;
00134             }
00135         case NearestToCand: {
00136             if (selVtxs_.empty() || (bestCand_ == 0)) return false;
00137             const reco::Vertex * which = 0;
00138             float dzmin = 9999.0; 
00139             for (vector<const reco::Vertex *>::const_iterator itv = selVtxs_.begin(), edv = selVtxs_.end(); itv != edv; ++itv) {
00140                 float dz = std::abs((*itv)->z() - bestCand_->vz());
00141                 if (dz < dzmin) { dzmin = dz; which = *itv; }
00142             }
00143             if (which == 0) return false; // actually it should not happen, but better safe than sorry
00144             auto_ptr<vector<reco::Vertex> > result(new vector<reco::Vertex>(1, *which));
00145             iEvent.put(result);
00146             return true;
00147             }
00148         case FromBeamSpot: {
00149             Handle<reco::BeamSpot> beamSpot;
00150             iEvent.getByLabel(beamSpot_, beamSpot);
00151             reco::Vertex bs(beamSpot->position(), beamSpot->rotatedCovariance3D(), 0, 0, 0);
00152             auto_ptr<vector<reco::Vertex> > result(new vector<reco::Vertex>(1, bs));
00153             iEvent.put(result);
00154             return true;
00155             }
00156         default:
00157             return false;
00158     }
00159 }
00160 
00161 #include "FWCore/Framework/interface/MakerMacros.h"
00162 DEFINE_FWK_MODULE(PATSingleVertexSelector);