CMS 3D CMS Logo

Classes | Public Member Functions | Private Member Functions | Private Attributes

CombinedSVComputer Class Reference

#include <CombinedSVComputer.h>

List of all members.

Classes

struct  IterationRange

Public Member Functions

 CombinedSVComputer (const edm::ParameterSet &params)
reco::TaggingVariableList operator() (const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const

Private Member Functions

IterationRange flipIterate (int size, bool vertex) const
double flipValue (double value, bool vertex) const
const
reco::TrackIPTagInfo::TrackIPData
threshTrack (const reco::TrackIPTagInfo &trackIPTagInfo, const reco::TrackIPTagInfo::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const

Private Attributes

double charmCut
double minTrackWeight
unsigned int pseudoMultiplicityMin
reco::V0Filter pseudoVertexV0Filter
reco::TrackIPTagInfo::SortCriteria sortCriterium
bool trackFlip
unsigned int trackMultiplicityMin
reco::TrackSelector trackNoDeltaRSelector
reco::V0Filter trackPairV0Filter
reco::TrackSelector trackPseudoSelector
reco::TrackSelector trackSelector
bool useTrackWeights
bool vertexFlip
bool vertexMassCorrection

Detailed Description

Definition at line 14 of file CombinedSVComputer.h.


Constructor & Destructor Documentation

CombinedSVComputer::CombinedSVComputer ( const edm::ParameterSet params)

Definition at line 51 of file CombinedSVComputer.cc.

                                                                    :
        trackFlip(params.getParameter<bool>("trackFlip")),
        vertexFlip(params.getParameter<bool>("vertexFlip")),
        charmCut(params.getParameter<double>("charmCut")),
        sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
        trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
        trackNoDeltaRSelector(dropDeltaR(params.getParameter<edm::ParameterSet>("trackSelection"))),
        trackPseudoSelector(params.getParameter<edm::ParameterSet>("trackPseudoSelection")),
        pseudoMultiplicityMin(params.getParameter<unsigned int>("pseudoMultiplicityMin")),
        trackMultiplicityMin(params.getParameter<unsigned int>("trackMultiplicityMin")),
        minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
        useTrackWeights(params.getParameter<bool>("useTrackWeights")),
        vertexMassCorrection(params.getParameter<bool>("correctVertexMass")),
        pseudoVertexV0Filter(params.getParameter<edm::ParameterSet>("pseudoVertexV0Filter")),
        trackPairV0Filter(params.getParameter<edm::ParameterSet>("trackPairV0Filter"))
{
}

Member Function Documentation

CombinedSVComputer::IterationRange CombinedSVComputer::flipIterate ( int  size,
bool  vertex 
) const [inline, private]

Definition at line 74 of file CombinedSVComputer.cc.

References CombinedSVComputer::IterationRange::begin, CombinedSVComputer::IterationRange::end, CombinedSVComputer::IterationRange::increment, findQualityFiles::size, trackFlip, and vertexFlip.

Referenced by operator()(), and threshTrack().

{
        IterationRange range;
        if (vertex ? vertexFlip : trackFlip) {
                range.begin = size - 1;
                range.end = -1;
                range.increment = -1;
        } else {
                range.begin = 0;
                range.end = size;
                range.increment = +1;
        }

        return range;
}
double CombinedSVComputer::flipValue ( double  value,
bool  vertex 
) const [inline, private]

Definition at line 69 of file CombinedSVComputer.cc.

References trackFlip, relativeConstraints::value, and vertexFlip.

Referenced by operator()().

{
        return (vertex ? vertexFlip : trackFlip) ? -value : value;
}
TaggingVariableList CombinedSVComputer::operator() ( const reco::TrackIPTagInfo ipInfo,
const reco::SecondaryVertexTagInfo svInfo 
) const

Definition at line 140 of file CombinedSVComputer.cc.

References edm::RefVector< C, T, F >::begin(), reco::TrackIPTagInfo::TrackIPData::closestToJetAxis, runTheMatrix::data, Geom::deltaR(), dir, reco::TrackIPTagInfo::TrackIPData::distanceToJetAxis, edm::RefVector< C, T, F >::end(), etaRel(), reco::SecondaryVertexTagInfo::flightDirection(), reco::SecondaryVertexTagInfo::flightDistance(), reco::btau::flightDistance2dSig, reco::btau::flightDistance2dVal, reco::btau::flightDistance3dSig, reco::btau::flightDistance3dVal, flipIterate(), flipValue(), reco::Vertex::hasRefittedTracks(), i, reco::TrackIPTagInfo::impactParameterData(), reco::TaggingVariableList::insert(), reco::TrackIPTagInfo::TrackIPData::ip2d, reco::TrackIPTagInfo::TrackIPData::ip3d, edm::Ref< C, T, F >::isNonnull(), j, reco::JTATagInfo::jet(), metsig::jet, reco::btau::jetEta, reco::btau::jetNSecondaryVertices, reco::btau::jetPt, PV3DBase< T, PVType, FrameType >::mag2(), minTrackWeight, reco::TrackBase::momentum(), reco::btag::Vertices::NoVertex, reco::SecondaryVertexTagInfo::nVertexTracks(), reco::SecondaryVertexTagInfo::nVertices(), convertSQLiteXML::ok, reco::TrackIPTagInfo::primaryVertex(), pseudoMultiplicityMin, reco::btag::Vertices::PseudoVertex, pseudoVertexV0Filter, range_for, reco::btag::Vertices::RecoVertex, reco::Vertex::refittedTrack(), reco::SecondaryVertexTagInfo::secondaryVertex(), reco::TrackIPTagInfo::selectedTracks(), Measurement1D::significance(), edm::RefVector< C, T, F >::size(), sortCriterium, reco::TrackIPTagInfo::sortedIndexes(), mathSSE::sqrt(), threshTrack(), reco::btau::trackDecayLenVal, reco::btau::trackDeltaR, reco::btau::trackEta, reco::btau::trackEtaRel, reco::btau::trackJetDistVal, reco::btau::trackMomentum, trackMultiplicityMin, trackPairV0Filter, reco::btau::trackPPar, reco::btau::trackPParRatio, trackPseudoSelector, reco::btau::trackPtRatio, reco::btau::trackPtRel, reco::JTATagInfo::tracks(), testEve_cfg::tracks, trackSelector, reco::btau::trackSip2dSig, reco::btau::trackSip2dSigAboveCharm, reco::btau::trackSip2dVal, reco::btau::trackSip3dSig, reco::btau::trackSip3dSigAboveCharm, reco::btau::trackSip3dVal, reco::btau::trackSumJetDeltaR, reco::btau::trackSumJetEtRatio, reco::SecondaryVertexTagInfo::trackWeight(), useTrackWeights, Measurement1D::value(), reco::btau::vertexCategory, reco::btau::vertexEnergyRatio, reco::btau::vertexJetDeltaR, reco::btau::vertexMass, vertexMassCorrection, reco::btau::vertexNTracks, reco::SecondaryVertexTagInfo::vertexTracks(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

{
        using namespace ROOT::Math;

        edm::RefToBase<Jet> jet = ipInfo.jet();
        math::XYZVector jetDir = jet->momentum().Unit();
        bool havePv = ipInfo.primaryVertex().isNonnull();
        GlobalPoint pv;
        if (havePv)
                pv = GlobalPoint(ipInfo.primaryVertex()->x(),
                                 ipInfo.primaryVertex()->y(),
                                 ipInfo.primaryVertex()->z());

        btag::Vertices::VertexType vtxType = btag::Vertices::NoVertex;

        TaggingVariableList vars; // = ipInfo.taggingVariables();

        vars.insert(btau::jetPt, jet->pt(), true);
        vars.insert(btau::jetEta, jet->eta(), true);

        if (ipInfo.tracks().size() < trackMultiplicityMin)
                return vars;

        TrackKinematics allKinematics;
        TrackKinematics vertexKinematics;

        int vtx = -1;
        IterationRange range = flipIterate(svInfo.nVertices(), true);
        range_for(i, range) {
                if (vtx < 0)
                        vtx = i;

                const Vertex &vertex = svInfo.secondaryVertex(i);
                bool hasRefittedTracks = vertex.hasRefittedTracks();
                TrackRefVector tracks = svInfo.vertexTracks(i);
                for(TrackRefVector::const_iterator track = tracks.begin();
                    track != tracks.end(); track++) {
                        double w = svInfo.trackWeight(i, *track);
                        if (w < minTrackWeight)
                                continue;
                        if (hasRefittedTracks) {
                                Track actualTrack =
                                                vertex.refittedTrack(*track);
                                vertexKinematics.add(actualTrack, w);
                                vars.insert(btau::trackEtaRel, etaRel(jetDir,
                                                actualTrack.momentum()), true);
                        } else {
                                vertexKinematics.add(**track, w);
                                vars.insert(btau::trackEtaRel, etaRel(jetDir,
                                                (*track)->momentum()), true);
                        }
                }
        }

        if (vtx >= 0) {
                vtxType = btag::Vertices::RecoVertex;

                vars.insert(btau::flightDistance2dVal,
                            flipValue(
                                svInfo.flightDistance(vtx, true).value(),
                                true),
                            true);
                vars.insert(btau::flightDistance2dSig,
                            flipValue(
                                svInfo.flightDistance(vtx, true).significance(),
                                true),
                            true);
                vars.insert(btau::flightDistance3dVal,
                            flipValue(
                                svInfo.flightDistance(vtx, false).value(),
                                true),
                            true);
                vars.insert(btau::flightDistance3dSig,
                            flipValue(
                                svInfo.flightDistance(vtx, false).significance(),
                                true),
                            true);
                vars.insert(btau::vertexJetDeltaR,
                            Geom::deltaR(svInfo.flightDirection(vtx), jetDir),true);
                vars.insert(btau::jetNSecondaryVertices, svInfo.nVertices(), true);
                vars.insert(btau::vertexNTracks, svInfo.nVertexTracks(), true);
        }

        std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
        const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
                                                ipInfo.impactParameterData();
        const edm::RefVector<TrackCollection> &tracks =
                                                ipInfo.selectedTracks();
        std::vector<TrackRef> pseudoVertexTracks;

        TrackRef trackPairV0Test[2];
        range = flipIterate(indices.size(), false);
        range_for(i, range) {
                std::size_t idx = indices[i];
                const TrackIPTagInfo::TrackIPData &data = ipData[idx];
                const TrackRef &trackRef = tracks[idx];
                const Track &track = *trackRef;

                // filter track

                if (!trackSelector(track, data, *jet, pv))
                        continue;

                // add track to kinematics for all tracks in jet

                allKinematics.add(track);

                // if no vertex was reconstructed, attempt pseudo vertex

                if (vtxType == btag::Vertices::NoVertex &&
                    trackPseudoSelector(track, data, *jet, pv)) {
                        pseudoVertexTracks.push_back(trackRef);
                        vertexKinematics.add(track);
                }

                // check against all other tracks for V0 track pairs

                trackPairV0Test[0] = tracks[idx];
                bool ok = true;
                range_for(j, range) {
                        if (i == j)
                                continue;

                        std::size_t pairIdx = indices[j];
                        const TrackIPTagInfo::TrackIPData &pairTrackData =
                                                        ipData[pairIdx];
                        const TrackRef &pairTrackRef = tracks[pairIdx];
                        const Track &pairTrack = *pairTrackRef;

                        if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
                                continue;

                        trackPairV0Test[1] = pairTrackRef;
                        if (!trackPairV0Filter(trackPairV0Test, 2)) {
                                ok = false;
                                break;
                        }
                }
                if (!ok)
                        continue;

                // add track variables

                math::XYZVector trackMom = track.momentum();
                double trackMag = std::sqrt(trackMom.Mag2());

                vars.insert(btau::trackSip3dVal,
                            flipValue(data.ip3d.value(), false), true);
                vars.insert(btau::trackSip3dSig,
                            flipValue(data.ip3d.significance(), false), true);
                vars.insert(btau::trackSip2dVal,
                            flipValue(data.ip2d.value(), false), true);
                vars.insert(btau::trackSip2dSig,
                            flipValue(data.ip2d.significance(), false), true);
                vars.insert(btau::trackJetDistVal, data.distanceToJetAxis.value(), true);
//              vars.insert(btau::trackJetDistSig, data.distanceToJetAxis.significance(), true);
//              vars.insert(btau::trackFirstTrackDist, data.distanceToFirstTrack, true);
//              vars.insert(btau::trackGhostTrackVal, data.distanceToGhostTrack.value(), true);
//              vars.insert(btau::trackGhostTrackSig, data.distanceToGhostTrack.significance(), true);
                vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToJetAxis - pv).mag() : -1.0, true);

                vars.insert(btau::trackMomentum, trackMag, true);
                vars.insert(btau::trackEta, trackMom.Eta(), true);

                vars.insert(btau::trackPtRel, VectorUtil::Perp(trackMom, jetDir), true);
                vars.insert(btau::trackPPar, jetDir.Dot(trackMom), true);
                vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
                vars.insert(btau::trackPtRatio, VectorUtil::Perp(trackMom, jetDir) / trackMag, true);
                vars.insert(btau::trackPParRatio, jetDir.Dot(trackMom) / trackMag, true);
        } 

        if (vtxType == btag::Vertices::NoVertex &&
            vertexKinematics.numberOfTracks() >= pseudoMultiplicityMin &&
            pseudoVertexV0Filter(pseudoVertexTracks)) {
                vtxType = btag::Vertices::PseudoVertex;
                for(std::vector<TrackRef>::const_iterator track =
                                                pseudoVertexTracks.begin();
                    track != pseudoVertexTracks.end(); ++track)
                        vars.insert(btau::trackEtaRel, etaRel(jetDir,
                                                (*track)->momentum()), true);
        }

        vars.insert(btau::vertexCategory, vtxType, true);
        vars.insert(btau::trackSumJetDeltaR,
                    VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
        vars.insert(btau::trackSumJetEtRatio,
                    allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
        vars.insert(btau::trackSip3dSigAboveCharm,
                    flipValue(
                        threshTrack(ipInfo, TrackIPTagInfo::IP3DSig, *jet, pv)
                                                        .ip3d.significance(),
                        false),
                    true);
        vars.insert(btau::trackSip2dSigAboveCharm,
                    flipValue(
                        threshTrack(ipInfo, TrackIPTagInfo::IP2DSig, *jet, pv)
                                                        .ip2d.significance(),
                        false),
                    true);

        if (vtxType != btag::Vertices::NoVertex) {
                math::XYZTLorentzVector allSum = useTrackWeights
                        ? allKinematics.weightedVectorSum()
                        : allKinematics.vectorSum();
                math::XYZTLorentzVector vertexSum = useTrackWeights
                        ? vertexKinematics.weightedVectorSum()
                        : vertexKinematics.vectorSum();

                if (vtxType != btag::Vertices::RecoVertex) {
                        vars.insert(btau::vertexNTracks,
                                    vertexKinematics.numberOfTracks(), true);
                        vars.insert(btau::vertexJetDeltaR,
                                    VectorUtil::DeltaR(vertexSum, jetDir), true);
                }

                double vertexMass = vertexSum.M();
                if (vtxType == btag::Vertices::RecoVertex &&
                    vertexMassCorrection) {
                        GlobalVector dir = svInfo.flightDirection(vtx);
                        double vertexPt2 =
                                math::XYZVector(dir.x(), dir.y(), dir.z()).
                                        Cross(vertexSum).Mag2() / dir.mag2();
                        vertexMass = std::sqrt(vertexMass * vertexMass +
                                               vertexPt2) + std::sqrt(vertexPt2);
                }
                vars.insert(btau::vertexMass, vertexMass, true);
                if (allKinematics.numberOfTracks())
                        vars.insert(btau::vertexEnergyRatio,
                                    vertexSum.E() / allSum.E(), true);
                else
                        vars.insert(btau::vertexEnergyRatio, 1, true);
        }

        vars.finalize();

        return vars;
}
const TrackIPTagInfo::TrackIPData & CombinedSVComputer::threshTrack ( const reco::TrackIPTagInfo trackIPTagInfo,
const reco::TrackIPTagInfo::SortCriteria  sort,
const reco::Jet jet,
const GlobalPoint pv 
) const [private]

Definition at line 92 of file CombinedSVComputer.cc.

References charmCut, runTheMatrix::data, flipIterate(), i, reco::TrackIPTagInfo::impactParameterData(), range_for, reco::TrackIPTagInfo::selectedTracks(), reco::TrackIPTagInfo::sortedIndexes(), trackNoDeltaRSelector, and testEve_cfg::tracks.

Referenced by operator()().

{
        const edm::RefVector<TrackCollection> &tracks =
                                        trackIPTagInfo.selectedTracks();
        const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
                                        trackIPTagInfo.impactParameterData();
        std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);

        IterationRange range = flipIterate(indices.size(), false);
        TrackKinematics kin;
        range_for(i, range) {
                std::size_t idx = indices[i];
                const TrackIPTagInfo::TrackIPData &data = ipData[idx];
                const Track &track = *tracks[idx];

                if (!trackNoDeltaRSelector(track, data, jet, pv))
                        continue;

                kin.add(track);
                if (kin.vectorSum().M() > charmCut) 
                        return data;
        }

        static const TrackIPTagInfo::TrackIPData dummy = {
                GlobalPoint(),
                GlobalPoint(),
                Measurement1D(-1.0, 1.0),
                Measurement1D(-1.0, 1.0),
                Measurement1D(-1.0, 1.0),
                Measurement1D(-1.0, 1.0),
                0.
        };
        return dummy;
}

Member Data Documentation

double CombinedSVComputer::charmCut [private]

Definition at line 36 of file CombinedSVComputer.h.

Referenced by threshTrack().

Definition at line 43 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 41 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 46 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 37 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 34 of file CombinedSVComputer.h.

Referenced by flipIterate(), and flipValue().

Definition at line 42 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 39 of file CombinedSVComputer.h.

Referenced by threshTrack().

Definition at line 47 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 40 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 38 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 44 of file CombinedSVComputer.h.

Referenced by operator()().

Definition at line 35 of file CombinedSVComputer.h.

Referenced by flipIterate(), and flipValue().

Definition at line 45 of file CombinedSVComputer.h.

Referenced by operator()().