![]() |
![]() |
#include <HLTVertexFilter.cc>
Public Member Functions | |
HLTVertexFilter (const edm::ParameterSet &config) | |
~HLTVertexFilter () | |
Private Member Functions | |
virtual bool | filter (edm::Event &event, const edm::EventSetup &setup) |
Private Attributes | |
edm::InputTag | m_inputTag |
double | m_maxChi2 |
double | m_maxD0 |
double | m_maxZ |
double | m_minNDoF |
unsigned int | m_minVertices |
Description: HLTFilter to accept events with at least a given number of vertices
Implementation: This class implements an HLTFilter to select events with at least a certain number of vertices matching some selection criteria.
Definition at line 36 of file HLTVertexFilter.cc.
HLTVertexFilter::HLTVertexFilter | ( | const edm::ParameterSet & | config | ) | [explicit] |
Definition at line 57 of file HLTVertexFilter.cc.
: m_inputTag(config.getParameter<edm::InputTag>("inputTag")), m_minNDoF(config.getParameter<double>("minNDoF")), m_maxChi2(config.getParameter<double>("maxChi2")), m_maxD0(config.getParameter<double>("maxD0")), m_maxZ(config.getParameter<double>("maxZ")), m_minVertices(config.getParameter<unsigned int>("minVertices")) { // register your products produces<trigger::TriggerFilterObjectWithRefs>(); }
HLTVertexFilter::~HLTVertexFilter | ( | ) |
Definition at line 70 of file HLTVertexFilter.cc.
{ }
bool HLTVertexFilter::filter | ( | edm::Event & | event, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
Implements HLTFilter.
Definition at line 81 of file HLTVertexFilter.cc.
References abs, reco::Vertex::chi2(), reco::Vertex::isFake(), reco::Vertex::isValid(), edm::HandleBase::isValid(), m_inputTag, m_maxChi2, m_maxD0, m_maxZ, m_minNDoF, m_minVertices, module(), n, reco::Vertex::ndof(), path(), funct::pow(), mathSSE::sqrt(), reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().
{ // The filter object std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterobject (new trigger::TriggerFilterObjectWithRefs(path(),module())); // get hold of collection of objects edm::Handle<reco::VertexCollection> vertices; event.getByLabel(m_inputTag, vertices); unsigned int n = 0; if (vertices.isValid()) { BOOST_FOREACH(const reco::Vertex & vertex, * vertices) { if (vertex.isValid() and not vertex.isFake() and (vertex.chi2() <= m_maxChi2) and (vertex.ndof() >= m_minNDoF) and (std::abs(vertex.z()) <= m_maxZ) and (std::abs(vertex.y()) <= m_maxD0) and (std::abs(vertex.x()) <= m_maxD0) and (std::sqrt(std::pow(vertex.x(),2)+std::pow(vertex.y(),2)) <= m_maxD0) ) ++n; } } // put filter object into the Event event.put(filterobject); // filter decision return (n >= m_minVertices); }
edm::InputTag HLTVertexFilter::m_inputTag [private] |
Definition at line 45 of file HLTVertexFilter.cc.
Referenced by filter().
double HLTVertexFilter::m_maxChi2 [private] |
Definition at line 47 of file HLTVertexFilter.cc.
Referenced by filter().
double HLTVertexFilter::m_maxD0 [private] |
Definition at line 48 of file HLTVertexFilter.cc.
Referenced by filter().
double HLTVertexFilter::m_maxZ [private] |
Definition at line 49 of file HLTVertexFilter.cc.
Referenced by filter().
double HLTVertexFilter::m_minNDoF [private] |
Definition at line 46 of file HLTVertexFilter.cc.
Referenced by filter().
unsigned int HLTVertexFilter::m_minVertices [private] |
Definition at line 50 of file HLTVertexFilter.cc.
Referenced by filter().