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 using namespace std;
00029
00030 modes_.push_back( parseMode(iConfig.getParameter<std::string>("mode")) );
00031 if (iConfig.exists("fallbacks")) {
00032 vector<string> modes = iConfig.getParameter<vector<string> >("fallbacks");
00033 for (vector<string>::const_iterator it = modes.begin(), ed = modes.end(); it != ed; ++it) {
00034 modes_.push_back( parseMode(*it) );
00035 }
00036 }
00037 if (hasMode_(First) || hasMode_(NearestToCand)) {
00038 vertices_ = iConfig.getParameter<edm::InputTag>("vertices");
00039 if (iConfig.existsAs<string>("vertexPreselection")) {
00040 string presel = iConfig.getParameter<string>("vertexPreselection");
00041 if (!presel.empty()) vtxPreselection_ = auto_ptr<VtxSel>(new VtxSel(presel));
00042 }
00043 }
00044 if (hasMode_(NearestToCand) || hasMode_(FromCand)) {
00045 candidates_ = iConfig.getParameter<vector<edm::InputTag> >("candidates");
00046 if (iConfig.existsAs<string>("candidatePreselection")) {
00047 string presel = iConfig.getParameter<string>("candidatePreselection");
00048 if (!presel.empty()) candPreselection_ = auto_ptr<CandSel>(new CandSel(presel));
00049 }
00050 }
00051 if (hasMode_(FromBeamSpot)) {
00052 beamSpot_ = iConfig.getParameter<edm::InputTag>("beamSpot");
00053 }
00054 produces<vector<reco::Vertex> >();
00055 }
00056
00057
00058 PATSingleVertexSelector::~PATSingleVertexSelector()
00059 {
00060 }
00061
00062 bool PATSingleVertexSelector::hasMode_(Mode mode) const {
00063 return (std::find(modes_.begin(), modes_.end(), mode) != modes_.end());
00064 }
00065
00066 bool
00067 PATSingleVertexSelector::filter(edm::Event & iEvent, const edm::EventSetup & iSetup) {
00068 using namespace edm;
00069 using namespace std;
00070
00071
00072 selVtxs_.clear(); bestCand_ = 0;
00073
00074
00075
00076 if (hasMode_(First) || hasMode_(NearestToCand)) {
00077 Handle<vector<reco::Vertex> > vertices;
00078 iEvent.getByLabel(vertices_, vertices);
00079 for (vector<reco::Vertex>::const_iterator itv = vertices->begin(), edv = vertices->end(); itv != edv; ++itv) {
00080 if ((vtxPreselection_.get() != 0) && !((*vtxPreselection_)(*itv)) ) continue;
00081 selVtxs_.push_back( &*itv );
00082 }
00083 }
00084
00085 if (hasMode_(NearestToCand) || hasMode_(FromCand)) {
00086 vector<pair<double, const reco::Candidate *> > cands;
00087 for (vector<edm::InputTag>::const_iterator itt = candidates_.begin(), edt = candidates_.end(); itt != edt; ++itt) {
00088 Handle<View<reco::Candidate> > theseCands;
00089 iEvent.getByLabel(*itt, theseCands);
00090 for (View<reco::Candidate>::const_iterator itc = theseCands->begin(), edc = theseCands->end(); itc != edc; ++itc) {
00091 if ((candPreselection_.get() != 0) && !((*candPreselection_)(*itc))) continue;
00092 cands.push_back( pair<double, const reco::Candidate *>(-itc->pt(), &*itc) );
00093 }
00094 }
00095 if (!cands.empty()) bestCand_ = cands.front().second;
00096 }
00097
00098
00099 for (std::vector<Mode>::const_iterator itm = modes_.begin(), endm = modes_.end(); itm != endm; ++itm) {
00100 if (filter_(*itm, iEvent, iSetup)) return true;
00101 }
00102 return false;
00103
00104
00105 selVtxs_.clear(); bestCand_ = 0;
00106 }
00107
00108 bool
00109 PATSingleVertexSelector::filter_(Mode mode, edm::Event & iEvent, const edm::EventSetup & iSetup) {
00110 using namespace edm;
00111 using namespace std;
00112 switch(mode) {
00113 case First: {
00114 if (selVtxs_.empty()) return false;
00115 auto_ptr<vector<reco::Vertex> > result(new vector<reco::Vertex>(1, *selVtxs_.front()));
00116 iEvent.put(result);
00117 return true;
00118 }
00119 case FromCand: {
00120 if (bestCand_ == 0) return false;
00121 reco::Vertex vtx;
00122 if (typeid(*bestCand_) == typeid(reco::VertexCompositeCandidate)) {
00123 vtx = reco::Vertex(bestCand_->vertex(), bestCand_->vertexCovariance(),
00124 bestCand_->vertexChi2(), bestCand_->vertexNdof(), bestCand_->numberOfDaughters() );
00125 } else {
00126 vtx = reco::Vertex(bestCand_->vertex(), reco::Vertex::Error(), 0, 0, 0);
00127 }
00128 auto_ptr<vector<reco::Vertex> > result(new vector<reco::Vertex>(1, vtx));
00129 iEvent.put(result);
00130 return true;
00131 }
00132 case NearestToCand: {
00133 if (selVtxs_.empty() || (bestCand_ == 0)) return false;
00134 const reco::Vertex * which = 0;
00135 float dzmin = 9999.0;
00136 for (vector<const reco::Vertex *>::const_iterator itv = selVtxs_.begin(), edv = selVtxs_.end(); itv != edv; ++itv) {
00137 float dz = abs((*itv)->z() - bestCand_->vz());
00138 if (dz < dzmin) { dzmin = dz; which = *itv; }
00139 }
00140 if (which == 0) return false;
00141 auto_ptr<vector<reco::Vertex> > result(new vector<reco::Vertex>(1, *which));
00142 iEvent.put(result);
00143 return true;
00144 }
00145 case FromBeamSpot: {
00146 Handle<reco::BeamSpot> beamSpot;
00147 iEvent.getByLabel(beamSpot_, beamSpot);
00148 reco::Vertex bs(beamSpot->position(), beamSpot->covariance3D(), 0, 0, 0);
00149 auto_ptr<vector<reco::Vertex> > result(new vector<reco::Vertex>(1, bs));
00150 iEvent.put(result);
00151 return true;
00152 }
00153 default:
00154 return false;
00155 }
00156 }
00157
00158 #include "FWCore/Framework/interface/MakerMacros.h"
00159 DEFINE_FWK_MODULE(PATSingleVertexSelector);