CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes
pat::PATSingleVertexSelector Class Reference

Produces a list containing a single vertex selected by some criteria. More...

#include "PhysicsTools/PatAlgos/plugins/PATSingleVertexSelector.h"

Inheritance diagram for pat::PATSingleVertexSelector:
edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

virtual bool filter (edm::Event &iEvent, const edm::EventSetup &iSetup)
 
 PATSingleVertexSelector (const edm::ParameterSet &iConfig)
 
 ~PATSingleVertexSelector ()
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
virtual ~EDFilter ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Types

typedef
StringCutObjectSelector
< reco::Candidate
CandSel
 
enum  Mode { First, NearestToCand, FromCand, FromBeamSpot }
 
typedef
StringCutObjectSelector
< reco::Vertex
VtxSel
 

Private Member Functions

std::auto_ptr< std::vector
< reco::Vertex > > 
filter_ (Mode mode, const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
bool hasMode_ (Mode mode) const
 

Static Private Member Functions

static Mode parseMode (const std::string &name)
 

Private Attributes

edm::InputTag beamSpot_
 
const reco::CandidatebestCand_
 
std::vector< edm::InputTagcandidates_
 
std::auto_ptr< CandSelcandPreselection_
 
bool doFilterEvents_
 
std::vector< Modemodes_
 
std::vector< const reco::Vertex * > selVtxs_
 
edm::InputTag vertices_
 
std::auto_ptr< VtxSelvtxPreselection_
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
typedef WorkerT< EDFilterWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDFilter
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Produces a list containing a single vertex selected by some criteria.

Author
Giovanni Petrucciani
Version
Id:
PATSingleVertexSelector.h,v 1.5 2011/06/15 11:47:25 friis Exp

Definition at line 27 of file PATSingleVertexSelector.h.

Member Typedef Documentation

Definition at line 39 of file PATSingleVertexSelector.h.

Definition at line 38 of file PATSingleVertexSelector.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

PATSingleVertexSelector::PATSingleVertexSelector ( const edm::ParameterSet iConfig)
explicit

Definition at line 27 of file PATSingleVertexSelector.cc.

References beamSpot_, candidates_, candPreselection_, doFilterEvents_, edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), First, FromBeamSpot, FromCand, edm::ParameterSet::getParameter(), hasMode_(), modes_, NearestToCand, parseMode(), vertices_, and vtxPreselection_.

28  : doFilterEvents_(false)
29 {
30  using namespace std;
31 
32  modes_.push_back( parseMode(iConfig.getParameter<std::string>("mode")) );
33  if (iConfig.exists("fallbacks")) {
34  vector<string> modes = iConfig.getParameter<vector<string> >("fallbacks");
35  for (vector<string>::const_iterator it = modes.begin(), ed = modes.end(); it != ed; ++it) {
36  modes_.push_back( parseMode(*it) );
37  }
38  }
40  vertices_ = iConfig.getParameter<edm::InputTag>("vertices");
41  if (iConfig.existsAs<string>("vertexPreselection")) {
42  string presel = iConfig.getParameter<string>("vertexPreselection");
43  if (!presel.empty()) vtxPreselection_ = auto_ptr<VtxSel>(new VtxSel(presel));
44  }
45  }
47  candidates_ = iConfig.getParameter<vector<edm::InputTag> >("candidates");
48  if (iConfig.existsAs<string>("candidatePreselection")) {
49  string presel = iConfig.getParameter<string>("candidatePreselection");
50  if (!presel.empty()) candPreselection_ = auto_ptr<CandSel>(new CandSel(presel));
51  }
52  }
53  if (hasMode_(FromBeamSpot)) {
54  beamSpot_ = iConfig.getParameter<edm::InputTag>("beamSpot");
55  }
56 
57  if ( iConfig.exists("filter") ) doFilterEvents_ = iConfig.getParameter<bool>("filter");
58 
59  produces<vector<reco::Vertex> >();
60 }
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:187
std::auto_ptr< VtxSel > vtxPreselection_
static Mode parseMode(const std::string &name)
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::auto_ptr< CandSel > candPreselection_
StringCutObjectSelector< reco::Vertex > VtxSel
StringCutObjectSelector< reco::Candidate > CandSel
std::vector< edm::InputTag > candidates_
PATSingleVertexSelector::~PATSingleVertexSelector ( )

Definition at line 63 of file PATSingleVertexSelector.cc.

64 {
65 }

Member Function Documentation

bool PATSingleVertexSelector::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDFilter.

Definition at line 72 of file PATSingleVertexSelector.cc.

References bestCand_, candidates_, candPreselection_, doFilterEvents_, filter_(), First, FromCand, edm::Event::getByLabel(), hasMode_(), modes_, NearestToCand, edm::Event::put(), query::result, selVtxs_, vertices_, and vtxPreselection_.

72  {
73  using namespace edm;
74  using namespace std;
75 
76  // Clear
77  selVtxs_.clear(); bestCand_ = 0;
78 
79  // Gather data from the Event
80  // -- vertex data --
82  Handle<vector<reco::Vertex> > vertices;
83  iEvent.getByLabel(vertices_, vertices);
84  for (vector<reco::Vertex>::const_iterator itv = vertices->begin(), edv = vertices->end(); itv != edv; ++itv) {
85  if ((vtxPreselection_.get() != 0) && !((*vtxPreselection_)(*itv)) ) continue;
86  selVtxs_.push_back( &*itv );
87  }
88  }
89  // -- candidate data --
91  vector<pair<double, const reco::Candidate *> > cands;
92  for (vector<edm::InputTag>::const_iterator itt = candidates_.begin(), edt = candidates_.end(); itt != edt; ++itt) {
93  Handle<View<reco::Candidate> > theseCands;
94  iEvent.getByLabel(*itt, theseCands);
95  for (View<reco::Candidate>::const_iterator itc = theseCands->begin(), edc = theseCands->end(); itc != edc; ++itc) {
96  if ((candPreselection_.get() != 0) && !((*candPreselection_)(*itc))) continue;
97  cands.push_back( pair<double, const reco::Candidate *>(-itc->pt(), &*itc) );
98  }
99  }
100  if (!cands.empty()) bestCand_ = cands.front().second;
101  }
102 
103  bool passes = false;
104  auto_ptr<vector<reco::Vertex> > result;
105  // Run main mode + possible fallback modes
106  for (std::vector<Mode>::const_iterator itm = modes_.begin(), endm = modes_.end(); itm != endm; ++itm) {
107  result = filter_(*itm, iEvent, iSetup);
108  // Check if we got any vertices. If so, take them.
109  if (result->size()) {
110  passes = true;
111  break;
112  }
113  }
114  iEvent.put(result);
115  // Check if we want to apply the EDFilter
116  if (doFilterEvents_)
117  return passes;
118  else return true;
119 }
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
std::auto_ptr< VtxSel > vtxPreselection_
std::auto_ptr< std::vector< reco::Vertex > > filter_(Mode mode, const edm::Event &iEvent, const edm::EventSetup &iSetup)
std::auto_ptr< CandSel > candPreselection_
std::vector< const reco::Vertex * > selVtxs_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
tuple result
Definition: query.py:137
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
std::vector< edm::InputTag > candidates_
const reco::Candidate * bestCand_
std::auto_ptr< std::vector< reco::Vertex > > PATSingleVertexSelector::filter_ ( Mode  mode,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)
private

Definition at line 122 of file PATSingleVertexSelector.cc.

References abs, SiPixelRawToDigiRegional_cfi::beamSpot, beamSpot_, bestCand_, First, FromBeamSpot, FromCand, edm::Event::getByLabel(), NearestToCand, reco::Candidate::numberOfDaughters(), query::result, selVtxs_, reco::Candidate::vertex(), reco::Candidate::vertexChi2(), reco::Candidate::vertexCovariance(), reco::Candidate::vertexNdof(), and reco::Candidate::vz().

Referenced by filter().

122  {
123  using namespace edm;
124  using namespace std;
125  std::auto_ptr<std::vector<reco::Vertex> > result(
126  new std::vector<reco::Vertex>());
127  switch(mode) {
128  case First: {
129  if (selVtxs_.empty()) return result;
130  result->push_back(*selVtxs_.front());
131  return result;
132  }
133  case FromCand: {
134  if (bestCand_ == 0) return result;
135  reco::Vertex vtx;
136  if (typeid(*bestCand_) == typeid(reco::VertexCompositeCandidate)) {
139  } else {
140  vtx = reco::Vertex(bestCand_->vertex(), reco::Vertex::Error(), 0, 0, 0);
141  }
142  result->push_back(vtx);
143  return result;
144  }
145  case NearestToCand: {
146  if (selVtxs_.empty() || (bestCand_ == 0)) return result;
147  const reco::Vertex * which = 0;
148  float dzmin = 9999.0;
149  for (vector<const reco::Vertex *>::const_iterator itv = selVtxs_.begin(), edv = selVtxs_.end(); itv != edv; ++itv) {
150  float dz = std::abs((*itv)->z() - bestCand_->vz());
151  if (dz < dzmin) { dzmin = dz; which = *itv; }
152  }
153  if (which != 0) // actually it should not happen, but better safe than sorry
154  result->push_back(*which);
155  return result;
156  }
157  case FromBeamSpot: {
159  iEvent.getByLabel(beamSpot_, beamSpot);
160  reco::Vertex bs(beamSpot->position(), beamSpot->covariance3D(), 0, 0, 0);
161  result->push_back(bs);
162  return result;
163  }
164  default:
165  // Return an empty vector signifying no vertices found.
166  return result;
167  }
168 }
#define abs(x)
Definition: mlp_lapack.h:159
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
virtual size_type numberOfDaughters() const =0
number of daughters
std::vector< const reco::Vertex * > selVtxs_
virtual double vertexChi2() const =0
chi-squares
tuple result
Definition: query.py:137
virtual const Point & vertex() const =0
vertex position
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual double vz() const =0
z coordinate of vertex position
virtual double vertexCovariance(int i, int j) const =0
(i, j)-th element of error matrix, i, j = 0, ... 2
const reco::Candidate * bestCand_
virtual double vertexNdof() const =0
bool PATSingleVertexSelector::hasMode_ ( Mode  mode) const
private

Definition at line 67 of file PATSingleVertexSelector.cc.

References spr::find(), alignBH_cfg::mode, and modes_.

Referenced by filter(), and PATSingleVertexSelector().

67  {
68  return (std::find(modes_.begin(), modes_.end(), mode) != modes_.end());
69 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
PATSingleVertexSelector::Mode PATSingleVertexSelector::parseMode ( const std::string &  name)
staticprivate

Definition at line 12 of file PATSingleVertexSelector.cc.

References edm::hlt::Exception, First, FromBeamSpot, FromCand, and NearestToCand.

Referenced by PATSingleVertexSelector().

12  {
13  if (mode == "firstVertex") {
14  return First;
15  } else if (mode == "nearestToCandidate") {
16  return NearestToCand;
17  } else if (mode == "fromCandidate") {
18  return FromCand;
19  } else if (mode == "beamSpot") {
20  return FromBeamSpot;
21  } else {
22  throw cms::Exception("Configuration") << "PATSingleVertexSelector: Mode '" << mode << "' not recognized or not supported.\n";
23  }
24 }

Member Data Documentation

edm::InputTag pat::PATSingleVertexSelector::beamSpot_
private

Definition at line 51 of file PATSingleVertexSelector.h.

Referenced by filter_(), and PATSingleVertexSelector().

const reco::Candidate* pat::PATSingleVertexSelector::bestCand_
private

Definition at line 54 of file PATSingleVertexSelector.h.

Referenced by filter(), and filter_().

std::vector<edm::InputTag> pat::PATSingleVertexSelector::candidates_
private

Definition at line 48 of file PATSingleVertexSelector.h.

Referenced by filter(), and PATSingleVertexSelector().

std::auto_ptr<CandSel> pat::PATSingleVertexSelector::candPreselection_
private

Definition at line 50 of file PATSingleVertexSelector.h.

Referenced by filter(), and PATSingleVertexSelector().

bool pat::PATSingleVertexSelector::doFilterEvents_
private

Definition at line 59 of file PATSingleVertexSelector.h.

Referenced by filter(), and PATSingleVertexSelector().

std::vector<Mode> pat::PATSingleVertexSelector::modes_
private

Definition at line 46 of file PATSingleVertexSelector.h.

Referenced by filter(), hasMode_(), and PATSingleVertexSelector().

std::vector<const reco::Vertex *> pat::PATSingleVertexSelector::selVtxs_
private

Definition at line 53 of file PATSingleVertexSelector.h.

Referenced by filter(), and filter_().

edm::InputTag pat::PATSingleVertexSelector::vertices_
private

Definition at line 47 of file PATSingleVertexSelector.h.

Referenced by filter(), and PATSingleVertexSelector().

std::auto_ptr<VtxSel > pat::PATSingleVertexSelector::vtxPreselection_
private

Definition at line 49 of file PATSingleVertexSelector.h.

Referenced by filter(), and PATSingleVertexSelector().