CMS 3D CMS Logo

Public Member Functions | Private Attributes

SimpleSecondaryVertexComputer Class Reference

#include <SimpleSecondaryVertexComputer.h>

Inheritance diagram for SimpleSecondaryVertexComputer:
JetTagComputer

List of all members.

Public Member Functions

float discriminator (const TagInfoHelper &tagInfos) const
 SimpleSecondaryVertexComputer (const edm::ParameterSet &parameters)

Private Attributes

unsigned int minTracks
bool unBoost
bool use2d
bool useSig

Detailed Description

Definition at line 16 of file SimpleSecondaryVertexComputer.h.


Constructor & Destructor Documentation

SimpleSecondaryVertexComputer::SimpleSecondaryVertexComputer ( const edm::ParameterSet parameters) [inline]

Definition at line 18 of file SimpleSecondaryVertexComputer.h.

References JetTagComputer::uses().

                                                                         :
                use2d(!parameters.getParameter<bool>("use3d")),
                useSig(parameters.getParameter<bool>("useSignificance")),
                unBoost(parameters.getParameter<bool>("unBoost")),
                minTracks(parameters.getParameter<unsigned int>("minTracks"))
        { uses("svTagInfos"); }

Member Function Documentation

float SimpleSecondaryVertexComputer::discriminator ( const TagInfoHelper tagInfos) const [inline, virtual]

Reimplemented from JetTagComputer.

Definition at line 25 of file SimpleSecondaryVertexComputer.h.

References reco::SecondaryVertexTagInfo::flightDistance(), JetTagComputer::TagInfoHelper::get(), info, funct::log(), minTracks, 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>();
                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;
        }

Member Data Documentation

Definition at line 65 of file SimpleSecondaryVertexComputer.h.

Referenced by discriminator().

Definition at line 64 of file SimpleSecondaryVertexComputer.h.

Referenced by discriminator().

Definition at line 62 of file SimpleSecondaryVertexComputer.h.

Referenced by discriminator().

Definition at line 63 of file SimpleSecondaryVertexComputer.h.

Referenced by discriminator().