CMS 3D CMS Logo

Classes | Functions

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoVertex/GhostTrackFitter/src/GhostTrackVertexFinder.cc File Reference

#include <algorithm>
#include <iterator>
#include <vector>
#include <map>
#include <set>
#include <Math/SVector.h>
#include <Math/SMatrix.h>
#include <Math/MatrixFunctions.h>
#include "DataFormats/GeometryVector/interface/GlobalPoint.h"
#include "DataFormats/GeometryVector/interface/GlobalVector.h"
#include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
#include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
#include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
#include "TrackingTools/GeomPropagators/interface/AnalyticalImpactPointExtrapolator.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
#include "TrackingTools/TransientTrack/interface/TransientTrackFromFTSFactory.h"
#include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
#include "RecoVertex/VertexPrimitives/interface/VertexFitter.h"
#include "RecoVertex/VertexPrimitives/interface/VertexState.h"
#include "RecoVertex/VertexPrimitives/interface/CachingVertex.h"
#include "RecoVertex/VertexPrimitives/interface/LinearizedTrackState.h"
#include "RecoVertex/VertexPrimitives/interface/ConvertError.h"
#include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
#include "RecoVertex/VertexTools/interface/LinearizedTrackStateFactory.h"
#include "RecoVertex/VertexTools/interface/VertexTrackFactory.h"
#include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
#include "RecoVertex/VertexTools/interface/GeometricAnnealing.h"
#include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
#include "RecoVertex/AdaptiveVertexFit/interface/AdaptiveVertexFitter.h"
#include "RecoVertex/GhostTrackFitter/interface/GhostTrack.h"
#include "RecoVertex/GhostTrackFitter/interface/GhostTrackState.h"
#include "RecoVertex/GhostTrackFitter/interface/GhostTrackFitter.h"
#include "RecoVertex/GhostTrackFitter/interface/GhostTrackPrediction.h"
#include "RecoVertex/GhostTrackFitter/interface/SequentialGhostTrackFitter.h"
#include "RecoVertex/GhostTrackFitter/interface/KalmanGhostTrackUpdater.h"
#include "RecoVertex/GhostTrackFitter/interface/GhostTrackVertexFinder.h"

Go to the source code of this file.

Classes

struct  reco::GhostTrackVertexFinder::FinderInfo

Functions

static bool covarianceUpdate (Matrix3S &cov, const Vector3 &residual, const Matrix3S &error, double &chi2, double theta, double phi)
static GhostTrackPrediction dummyPrediction (const Vertex &primaryVertex, const Track &ghostTrack)
static double fitChi2 (const CachingVertex< 5 > &vtx)
static void mergeTrackHelper (const std::vector< RefCountedVertexTrack > &tracks, std::vector< RefCountedVertexTrack > &newTracks, const VertexState &state, const VtxTrackIs &ghostTrackFinder, RefCountedVertexTrack &ghostTrack, const VertexTrackFactory< 5 > &factory)
static RefCountedVertexTrack relinearizeTrack (const RefCountedVertexTrack &track, const VertexState &state, const VertexTrackFactory< 5 > factory)
static std::vector
< RefCountedVertexTrack > 
relinearizeTracks (const std::vector< RefCountedVertexTrack > &tracks, const VertexState &state)
static VertexState stateMean (const VertexState &v1, const VertexState &v2)
static double trackVertexCompat (const CachingVertex< 5 > &vtx, const RefCountedVertexTrack &vertexTrack)
static TransientTrack transientGhostTrack (const GhostTrackPrediction &pred, const MagneticField *field)
static CachingVertex< 5 > vertexAtState (const TransientTrack &ghostTrack, const GhostTrackPrediction &pred, const GhostTrackState &state)
static double vtxErrorLong (const GlobalError &error, const GlobalVector &dir)
static GlobalPoint vtxMean (const GlobalPoint &p1, const GlobalError &e1, const GlobalPoint &p2, const GlobalError &e2)

Function Documentation

static bool covarianceUpdate ( Matrix3S &  cov,
const Vector3 &  residual,
const Matrix3S &  error,
double &  chi2,
double  theta,
double  phi 
) [static]

Definition at line 139 of file GhostTrackVertexFinder.cc.

References funct::cos(), and funct::sin().

Referenced by vertexAtState().

{
        using namespace ROOT::Math;

        Matrix23 jacobian;
        jacobian(0, 0) = std::cos(phi) * std::cos(theta);
        jacobian(0, 1) = std::sin(phi) * std::cos(theta);
        jacobian(0, 2) = -std::sin(theta);
        jacobian(1, 0) = -std::sin(phi);
        jacobian(1, 1) = std::cos(phi);

        Matrix2S measErr = Similarity(jacobian, error);
        Matrix2S combErr = Similarity(jacobian, cov) + measErr;
        if (!measErr.Invert() || !combErr.Invert())
                return false;

        cov -= SimilarityT(jacobian * cov, combErr);
        chi2 += Similarity(jacobian * residual, measErr);

        return true;
}
static GhostTrackPrediction dummyPrediction ( const Vertex primaryVertex,
const Track ghostTrack 
) [static]
static double fitChi2 ( const CachingVertex< 5 > &  vtx) [static]
static void mergeTrackHelper ( const std::vector< RefCountedVertexTrack > &  tracks,
std::vector< RefCountedVertexTrack > &  newTracks,
const VertexState state,
const VtxTrackIs &  ghostTrackFinder,
RefCountedVertexTrack &  ghostTrack,
const VertexTrackFactory< 5 > &  factory 
) [static]

Definition at line 777 of file GhostTrackVertexFinder.cc.

References gt, and relinearizeTrack().

Referenced by reco::GhostTrackVertexFinder::mergeVertices().

{
        for(std::vector<RefCountedVertexTrack>::const_iterator iter =
                        tracks.begin(); iter != tracks.end(); ++iter) {
                bool gt = ghostTrackFinder(*iter);
                if (gt && ghostTrack)
                        continue;

                RefCountedVertexTrack track =
                                relinearizeTrack(*iter, state, factory);

                if (gt)
                        ghostTrack = *iter;
                else
                        newTracks.push_back(*iter);
        }
}
static RefCountedVertexTrack relinearizeTrack ( const RefCountedVertexTrack &  track,
const VertexState state,
const VertexTrackFactory< 5 >  factory 
) [static]

Definition at line 211 of file GhostTrackVertexFinder.cc.

References VertexState::position(), and VertexTrackFactory< N >::vertexTrack().

Referenced by mergeTrackHelper(), reco::GhostTrackVertexFinder::reassignTracks(), and relinearizeTracks().

{
        RefCountedLinearizedTrackState linTrack = track->linearizedTrack();
        linTrack = linTrack->stateWithNewLinearizationPoint(state.position());
        return factory.vertexTrack(linTrack, state);
}
static std::vector<RefCountedVertexTrack> relinearizeTracks ( const std::vector< RefCountedVertexTrack > &  tracks,
const VertexState state 
) [static]

Definition at line 221 of file GhostTrackVertexFinder.cc.

References relinearizeTrack().

Referenced by reco::GhostTrackVertexFinder::reassignTracks().

{
        VertexTrackFactory<5> vertexTrackFactory;

        std::vector<RefCountedVertexTrack> finalTracks;
        finalTracks.reserve(tracks.size());

        for(std::vector<RefCountedVertexTrack>::const_iterator iter =
                                tracks.begin(); iter != tracks.end(); ++iter)
                finalTracks.push_back(
                        relinearizeTrack(*iter, state, vertexTrackFactory));

        return finalTracks;
}
static VertexState stateMean ( const VertexState v1,
const VertexState v2 
) [static]

Definition at line 132 of file GhostTrackVertexFinder.cc.

References VertexState::error(), VertexState::position(), and vtxMean().

Referenced by reco::GhostTrackVertexFinder::mergeVertices().

{
        return VertexState(vtxMean(v1.position(), v1.error(),
                                   v2.position(), v2.error()),
                           v1.error() + v2.error());
}
static double trackVertexCompat ( const CachingVertex< 5 > &  vtx,
const RefCountedVertexTrack &  vertexTrack 
) [static]

Definition at line 250 of file GhostTrackVertexFinder.cc.

References conv, dir, CachingVertex< N >::error(), error, GlobalErrorBase< T, ErrorWeightType >::matrix_new(), point, CachingVertex< N >::position(), and reco::TransientTrack::track().

Referenced by reco::GhostTrackVertexFinder::reassignTracks().

{
        using namespace ROOT::Math;

        TransientTrack track = vertexTrack->linearizedTrack()->track();
        GlobalPoint point = vtx.position();
        AnalyticalImpactPointExtrapolator extrap(track.field());
        TrajectoryStateOnSurface tsos =
                        extrap.extrapolate(track.impactPointState(), point);

        if (!tsos.isValid())
                return 1.0e6;

        GlobalPoint point1 = vtx.position();
        GlobalPoint point2 = tsos.globalPosition();
        Vector3 dir = conv(point2 - point1);
        Matrix3S error = vtx.error().matrix_new() +
                         tsos.cartesianError().matrix().Sub<Matrix3S>(0, 0);
        if (!error.Invert())
                return 1.0e6;

        return ROOT::Math::Similarity(error, dir);
}
static TransientTrack transientGhostTrack ( const GhostTrackPrediction pred,
const MagneticField field 
) [static]
static CachingVertex<5> vertexAtState ( const TransientTrack ghostTrack,
const GhostTrackPrediction pred,
const GhostTrackState state 
) [static]

Definition at line 163 of file GhostTrackVertexFinder.cc.

References reco::GhostTrackState::cartesianError(), conv, covarianceUpdate(), error, reco::GhostTrackState::globalPosition(), reco::GhostTrackState::isValid(), reco::GhostTrackState::lambda(), LinearizedTrackStateFactory::linearizedTrackState(), GlobalErrorBase< T, ErrorWeightType >::matrix_new(), point, reco::GhostTrackPrediction::position(), reco::GhostTrackPrediction::positionError(), reco::GhostTrackState::track(), VertexTrackFactory< N >::vertexTrack(), and vtxMean().

Referenced by reco::GhostTrackVertexFinder::initialVertices(), reco::GhostTrackVertexFinder::reassignTracks(), and reco::GhostTrackVertexFinder::refitGhostTrack().

{
        LinearizedTrackStateFactory linTrackFactory;
        VertexTrackFactory<5> vertexTrackFactory;

        if (!state.isValid())
                return CachingVertex<5>();

        GlobalPoint pca1 = pred.position(state.lambda());
        GlobalError err1 = pred.positionError(state.lambda());

        GlobalPoint pca2 = state.globalPosition();
        GlobalError err2 = state.cartesianError();

        GlobalPoint point = vtxMean(pca1, err1, pca2, err2);

        TransientTrack recTrack = state.track();

        RefCountedLinearizedTrackState linState[2] = {
                linTrackFactory.linearizedTrackState(point, ghostTrack),
                linTrackFactory.linearizedTrackState(point, recTrack)
        };
        if ( !linState[0]->isValid() || !linState[1]->isValid() )
          return CachingVertex<5>();

        Matrix3S cov = SMatrixIdentity();
        cov *= 10000;
        double chi2 = 0.;
        if (!covarianceUpdate(cov, conv(pca1 - point), err1.matrix_new(), chi2,
                               linState[0]->predictedStateParameters()[1],
                               linState[0]->predictedStateParameters()[2]) ||
            !covarianceUpdate(cov, conv(pca2 - point), err2.matrix_new(), chi2,
                               linState[1]->predictedStateParameters()[1],
                               linState[1]->predictedStateParameters()[2]))
                return CachingVertex<5>();

        GlobalError error(cov);
        VertexState vtxState(point, error);

        std::vector<RefCountedVertexTrack> linTracks(2);
        linTracks[0] = vertexTrackFactory.vertexTrack(linState[0], vtxState);
        linTracks[1] = vertexTrackFactory.vertexTrack(linState[1], vtxState);

        return CachingVertex<5>(point, error, linTracks, chi2);
}
static double vtxErrorLong ( const GlobalError error,
const GlobalVector dir 
) [static]

Definition at line 114 of file GhostTrackVertexFinder.cc.

References conv, and GlobalErrorBase< T, ErrorWeightType >::matrix_new().

Referenced by vtxMean().

{
        return ROOT::Math::Similarity(conv(dir), error.matrix_new());
}
static GlobalPoint vtxMean ( const GlobalPoint p1,
const GlobalError e1,
const GlobalPoint p2,
const GlobalError e2 
) [static]

Definition at line 119 of file GhostTrackVertexFinder.cc.

References diffTreeTool::diff, p1, vtxErrorLong(), and CommonMethods::weight().

Referenced by GlobalTrackingRegion::checkRZ(), RectangularEtaPhiTrackingRegion::checkRZOld(), stateMean(), and vertexAtState().

{
        GlobalVector diff = p2 - p1;

        double err1 = vtxErrorLong(e1, diff);
        double err2 = vtxErrorLong(e2, diff);

        double weight = err1 / (err1 + err2);

        return p1 + weight * diff;
}