CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATSingleVertexSelector.cc
Go to the documentation of this file.
5 
7 
8 #include <algorithm>
9 
11 
12 
14 PATSingleVertexSelector::parseMode(const std::string &mode) const {
15  if (mode == "firstVertex") {
16  return First;
17  } else if (mode == "nearestToCandidate") {
18  return NearestToCand;
19  } else if (mode == "fromCandidate") {
20  return FromCand;
21  } else if (mode == "beamSpot") {
22  return FromBeamSpot;
23  } else {
24  throw cms::Exception("Configuration") << "PATSingleVertexSelector: Mode '" << mode << "' not recognized or not supported.\n";
25  }
26 }
27 
28 
30  vtxPreselection_( iConfig.existsAs<std::string>("vertexPreselection") ? iConfig.getParameter<std::string>("vertexPreselection") : std::string(" 1 == 1 ") ),
31  candPreselection_( iConfig.existsAs<std::string>("candidatePreselection") ? iConfig.getParameter<std::string>("candidatePreselection") : std::string(" 1 == 1 ") ),
32  doFilterEvents_(false)
33 {
34  using namespace std;
35 
36  modes_.push_back( parseMode(iConfig.getParameter<std::string>("mode")) );
37  if (iConfig.exists("fallbacks")) {
38  vector<string> modes = iConfig.getParameter<vector<string> >("fallbacks");
39  for (vector<string>::const_iterator it = modes.begin(), ed = modes.end(); it != ed; ++it) {
40  modes_.push_back( parseMode(*it) );
41  }
42  }
44  verticesToken_ = consumes<vector<reco::Vertex> >(iConfig.getParameter<edm::InputTag>("vertices"));
45  }
47  candidatesToken_ = edm::vector_transform(iConfig.getParameter<vector<edm::InputTag> >("candidates"), [this](edm::InputTag const & tag){return consumes<edm::View<reco::Candidate> >(tag);});
48  }
49  if (hasMode_(FromBeamSpot)) {
50  beamSpotToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"));
51  }
52 
53  if ( iConfig.exists("filter") ) doFilterEvents_ = iConfig.getParameter<bool>("filter");
54 
55  produces<vector<reco::Vertex> >();
56 }
57 
58 
60 {
61 }
62 
64  return (std::find(modes_.begin(), modes_.end(), mode) != modes_.end());
65 }
66 
67 bool
69  using namespace edm;
70  using namespace std;
71 
72  // Clear
74 
75  // Gather data from the Event
76  // -- vertex data --
79  iEvent.getByToken(verticesToken_, vertices);
80  for (vector<reco::Vertex>::const_iterator itv = vertices->begin(), edv = vertices->end(); itv != edv; ++itv) {
81  if ( !(vtxPreselection_(*itv)) ) continue;
82  selVtxs_.push_back( reco::VertexRef(vertices,std::distance(vertices->begin(),itv) ) );
83  }
84  }
85  // -- candidate data --
87  vector<pair<double, reco::CandidatePtr> > cands;
88  for (vector<edm::EDGetTokenT<edm::View<reco::Candidate> > >::const_iterator itt = candidatesToken_.begin(), edt = candidatesToken_.end(); itt != edt; ++itt) {
89  Handle<View<reco::Candidate> > theseCands;
90  iEvent.getByToken(*itt, theseCands);
91  for (View<reco::Candidate>::const_iterator itc = theseCands->begin(), edc = theseCands->end(); itc != edc; ++itc) {
92  if ( !(candPreselection_(*itc)) ) continue;
93  cands.push_back( pair<double, reco::CandidatePtr>(-itc->pt(), reco::CandidatePtr(theseCands,std::distance(theseCands->begin(),itc) ) ) );
94  }
95  }
96  if (!cands.empty()) bestCand_ = cands.front().second;
97  }
98 
99  bool passes = false;
100  auto_ptr<vector<reco::Vertex> > result;
101  // Run main mode + possible fallback modes
102  for (std::vector<Mode>::const_iterator itm = modes_.begin(), endm = modes_.end(); itm != endm; ++itm) {
103  result = filter_(*itm, iEvent, iSetup);
104  // Check if we got any vertices. If so, take them.
105  if (result->size()) {
106  passes = true;
107  break;
108  }
109  }
110  iEvent.put(result);
111  // Check if we want to apply the EDFilter
112  if (doFilterEvents_)
113  return passes;
114  else return true;
115 }
116 
117 std::auto_ptr<std::vector<reco::Vertex> >
119  using namespace edm;
120  using namespace std;
121  std::auto_ptr<std::vector<reco::Vertex> > result(new std::vector<reco::Vertex>());
122  switch(mode) {
123  case First: {
124  if (selVtxs_.empty()) return result;
125  result->push_back(*selVtxs_.front());
126  return result;
127  }
128  case FromCand: {
129  if (bestCand_.isNull()) return result;
130  reco::Vertex vtx;
131  if (typeid(*bestCand_) == typeid(reco::VertexCompositeCandidate)) {
132  vtx = reco::Vertex(bestCand_->vertex(), bestCand_->vertexCovariance(),
133  bestCand_->vertexChi2(), bestCand_->vertexNdof(), bestCand_->numberOfDaughters() );
134  } else {
135  vtx = reco::Vertex(bestCand_->vertex(), reco::Vertex::Error(), 0, 0, 0);
136  }
137  result->push_back(vtx);
138  return result;
139  }
140  case NearestToCand: {
141  if (selVtxs_.empty() || (bestCand_.isNull())) return result;
143  float dzmin = 9999.0;
144  for (auto itv = selVtxs_.begin(), edv = selVtxs_.end(); itv != edv; ++itv) {
145  float dz = std::abs((*itv)->z() - bestCand_->vz());
146  if (dz < dzmin) { dzmin = dz; which = *itv; }
147  }
148  if (which.isNonnull()) // actually it should not happen, but better safe than sorry
149  result->push_back(*which);
150  return result;
151  }
152  case FromBeamSpot: {
154  iEvent.getByToken(beamSpotToken_, beamSpot);
155  reco::Vertex bs(beamSpot->position(), beamSpot->covariance3D(), 0, 0, 0);
156  result->push_back(bs);
157  return result;
158  }
159  default:
160  // Return an empty vector signifying no vertices found.
161  return result;
162  }
163 }
164 
T getParameter(std::string const &) const
Produces a list containing a single vertex selected by some criteria.
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
def which
Definition: eostools.py:335
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::auto_ptr< std::vector< reco::Vertex > > filter_(Mode mode, const edm::Event &iEvent, const edm::EventSetup &iSetup)
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< std::vector< reco::Vertex > > verticesToken_
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< reco::VertexRef > selVtxs_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
tuple result
Definition: query.py:137
Mode parseMode(const std::string &name) const
bool isNull() const
Checks for null.
Definition: Ptr.h:165
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< edm::EDGetTokenT< edm::View< reco::Candidate > > > candidatesToken_
virtual bool filter(edm::Event &iEvent, const edm::EventSetup &iSetup) override
PATSingleVertexSelector(const edm::ParameterSet &iConfig)
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
volatile std::atomic< bool > shutdown_flag false