#include <SimpleSecondaryVertexComputer.h>
Public Member Functions | |
float | discriminator (const TagInfoHelper &tagInfos) const |
SimpleSecondaryVertexComputer (const edm::ParameterSet ¶meters) | |
Private Attributes | |
unsigned int | minTracks |
unsigned int | minVertices_ |
bool | unBoost |
bool | use2d |
bool | useSig |
Definition at line 16 of file SimpleSecondaryVertexComputer.h.
SimpleSecondaryVertexComputer::SimpleSecondaryVertexComputer | ( | const edm::ParameterSet & | parameters | ) | [inline] |
Definition at line 18 of file SimpleSecondaryVertexComputer.h.
References edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), minVertices_, and JetTagComputer::uses().
: use2d(!parameters.getParameter<bool>("use3d")), useSig(parameters.getParameter<bool>("useSignificance")), unBoost(parameters.getParameter<bool>("unBoost")), minTracks(parameters.getParameter<unsigned int>("minTracks")), minVertices_(1) { uses("svTagInfos"); minVertices_ = parameters.existsAs<unsigned int>("minVertices") ? parameters.getParameter<unsigned int>("minVertices") : 1 ; }
float SimpleSecondaryVertexComputer::discriminator | ( | const TagInfoHelper & | tagInfos | ) | const [inline, virtual] |
Reimplemented from JetTagComputer.
Definition at line 29 of file SimpleSecondaryVertexComputer.h.
References reco::SecondaryVertexTagInfo::flightDistance(), JetTagComputer::TagInfoHelper::get(), info, create_public_lumi_plots::log, minTracks, minVertices_, reco::SecondaryVertexTagInfo::nVertexTracks(), reco::SecondaryVertexTagInfo::nVertices(), reco::SecondaryVertexTagInfo::secondaryVertex(), Measurement1D::significance(), unBoost, use2d, useSig, Measurement1D::value(), relativeConstraints::value, and reco::TrackKinematics::vectorSum().
{ const reco::SecondaryVertexTagInfo &info = tagInfos.get<reco::SecondaryVertexTagInfo>(); if(info.nVertices() < minVertices_) return -1; unsigned int idx = 0; while(idx < info.nVertices()) { if (info.nVertexTracks(idx) >= minTracks) break; idx++; } if (idx >= info.nVertices()) return -1.0; double gamma; if (unBoost) { reco::TrackKinematics kinematics( info.secondaryVertex(idx)); gamma = kinematics.vectorSum().Gamma(); } else gamma = 1.0; double value; if (useSig) value = info.flightDistance(idx, use2d).significance(); else value = info.flightDistance(idx, use2d).value(); value /= gamma; if (useSig) value = (value > 0) ? +std::log(1 + value) : -std::log(1 - value); return value; }
unsigned int SimpleSecondaryVertexComputer::minTracks [private] |
Definition at line 70 of file SimpleSecondaryVertexComputer.h.
Referenced by discriminator().
unsigned int SimpleSecondaryVertexComputer::minVertices_ [private] |
Definition at line 71 of file SimpleSecondaryVertexComputer.h.
Referenced by discriminator(), and SimpleSecondaryVertexComputer().
bool SimpleSecondaryVertexComputer::unBoost [private] |
Definition at line 69 of file SimpleSecondaryVertexComputer.h.
Referenced by discriminator().
bool SimpleSecondaryVertexComputer::use2d [private] |
Definition at line 67 of file SimpleSecondaryVertexComputer.h.
Referenced by discriminator().
bool SimpleSecondaryVertexComputer::useSig [private] |
Definition at line 68 of file SimpleSecondaryVertexComputer.h.
Referenced by discriminator().