CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/RecoBTag/SecondaryVertex/plugins/SecondaryVertexProducer.cc

Go to the documentation of this file.
00001 #include <functional>
00002 #include <algorithm>
00003 #include <iterator>
00004 #include <cstddef>
00005 #include <string>
00006 #include <vector>
00007 #include <map>
00008 #include <set>
00009 
00010 #include <boost/iterator/transform_iterator.hpp>
00011 
00012 #include "FWCore/Framework/interface/EDProducer.h"
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Framework/interface/MakerMacros.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/Utilities/interface/InputTag.h"
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 #include "FWCore/Utilities/interface/Exception.h"
00020 
00021 #include "DataFormats/TrackReco/interface/Track.h"
00022 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00023 #include "DataFormats/VertexReco/interface/Vertex.h"
00024 
00025 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
00026 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
00027 #include "DataFormats/BTauReco/interface/SecondaryVertexTagInfo.h"
00028 
00029 #include "RecoVertex/VertexPrimitives/interface/VertexException.h"
00030 #include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
00031 #include "RecoVertex/ConfigurableVertexReco/interface/ConfigurableVertexReconstructor.h"
00032 #include "RecoVertex/GhostTrackFitter/interface/GhostTrackVertexFinder.h"
00033 #include "RecoVertex/GhostTrackFitter/interface/GhostTrackPrediction.h"
00034 #include "RecoVertex/GhostTrackFitter/interface/GhostTrackState.h"
00035 #include "RecoVertex/GhostTrackFitter/interface/GhostTrack.h"
00036 
00037 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00038 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00039 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00040 
00041 #include "RecoBTag/SecondaryVertex/interface/TrackSelector.h"
00042 #include "RecoBTag/SecondaryVertex/interface/TrackSorting.h"
00043 #include "RecoBTag/SecondaryVertex/interface/SecondaryVertex.h"
00044 #include "RecoBTag/SecondaryVertex/interface/VertexFilter.h"
00045 #include "RecoBTag/SecondaryVertex/interface/VertexSorting.h"
00046 
00047 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
00048 
00049 using namespace reco;
00050 
00051 namespace {
00052         template<typename T>
00053         struct RefToBaseLess : public std::binary_function<edm::RefToBase<T>,
00054                                                            edm::RefToBase<T>,
00055                                                            bool> {
00056                 inline bool operator()(const edm::RefToBase<T> &r1,
00057                                        const edm::RefToBase<T> &r2) const
00058                 {
00059                         return r1.id() < r2.id() ||
00060                                (r1.id() == r2.id() && r1.key() < r2.key());
00061                 }
00062         };
00063 }
00064 
00065 class SecondaryVertexProducer : public edm::EDProducer {
00066     public:
00067         explicit SecondaryVertexProducer(const edm::ParameterSet &params);
00068         ~SecondaryVertexProducer();
00069 
00070         virtual void produce(edm::Event &event, const edm::EventSetup &es);
00071 
00072     private:
00073         enum ConstraintType {
00074                 CONSTRAINT_NONE = 0,
00075                 CONSTRAINT_BEAMSPOT,
00076                 CONSTRAINT_PV_BEAMSPOT_SIZE,
00077                 CONSTRAINT_PV_BS_Z_ERRORS_SCALED,
00078                 CONSTRAINT_PV_ERROR_SCALED,
00079                 CONSTRAINT_PV_PRIMARIES_IN_FIT
00080         };
00081 
00082         static ConstraintType getConstraintType(const std::string &name);
00083 
00084         const edm::InputTag             trackIPTagInfoLabel;
00085         edm::InputTag beamSpotTag;
00086         TrackIPTagInfo::SortCriteria    sortCriterium;
00087         TrackSelector                   trackSelector;
00088         ConstraintType                  constraint;
00089         double                          constraintScaling;
00090         edm::ParameterSet               vtxRecoPSet;
00091         bool                            useGhostTrack;
00092         bool                            withPVError;
00093         double                          minTrackWeight;
00094         VertexFilter                    vertexFilter;
00095         VertexSorting                   vertexSorting;
00096         bool                            useExternalSV;  
00097         double                          extSVDeltaRToJet;      
00098         edm::InputTag                   extSVCollection; 
00099 };
00100 
00101 SecondaryVertexProducer::ConstraintType
00102 SecondaryVertexProducer::getConstraintType(const std::string &name)
00103 {
00104         if (name == "None")
00105                 return CONSTRAINT_NONE;
00106         else if (name == "BeamSpot")
00107                 return CONSTRAINT_BEAMSPOT;
00108         else if (name == "BeamSpot+PVPosition")
00109                 return CONSTRAINT_PV_BEAMSPOT_SIZE;
00110         else if (name == "BeamSpotZ+PVErrorScaledXY")
00111                 return CONSTRAINT_PV_BS_Z_ERRORS_SCALED;
00112         else if (name == "PVErrorScaled")
00113                 return CONSTRAINT_PV_ERROR_SCALED;
00114         else if (name == "BeamSpot+PVTracksInFit")
00115                 return CONSTRAINT_PV_PRIMARIES_IN_FIT;
00116         else
00117                 throw cms::Exception("InvalidArgument")
00118                         << "SecondaryVertexProducer: ``constraint'' parameter "
00119                            "value \"" << name << "\" not understood."
00120                         << std::endl;
00121 }
00122 
00123 static GhostTrackVertexFinder::FitType
00124 getGhostTrackFitType(const std::string &name)
00125 {
00126         if (name == "AlwaysWithGhostTrack")
00127                 return GhostTrackVertexFinder::kAlwaysWithGhostTrack;
00128         else if (name == "SingleTracksWithGhostTrack")
00129                 return GhostTrackVertexFinder::kSingleTracksWithGhostTrack;
00130         else if (name == "RefitGhostTrackWithVertices")
00131                 return GhostTrackVertexFinder::kRefitGhostTrackWithVertices;
00132         else
00133                 throw cms::Exception("InvalidArgument")
00134                         << "SecondaryVertexProducer: ``fitType'' "
00135                            "parameter value \"" << name << "\" for "
00136                            "GhostTrackVertexFinder settings not "
00137                            "understood." << std::endl;
00138 }
00139 
00140 SecondaryVertexProducer::SecondaryVertexProducer(
00141                                         const edm::ParameterSet &params) :
00142         trackIPTagInfoLabel(params.getParameter<edm::InputTag>("trackIPTagInfos")),
00143         sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
00144         trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
00145         constraint(getConstraintType(params.getParameter<std::string>("constraint"))),
00146         constraintScaling(1.0),
00147         vtxRecoPSet(params.getParameter<edm::ParameterSet>("vertexReco")),
00148         useGhostTrack(vtxRecoPSet.getParameter<std::string>("finder") == "gtvr"),
00149         withPVError(params.getParameter<bool>("usePVError")),
00150         minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
00151         vertexFilter(params.getParameter<edm::ParameterSet>("vertexCuts")),
00152         vertexSorting(params.getParameter<edm::ParameterSet>("vertexSelection"))
00153 {
00154         if (constraint == CONSTRAINT_PV_ERROR_SCALED ||
00155             constraint == CONSTRAINT_PV_BS_Z_ERRORS_SCALED)
00156                 constraintScaling = params.getParameter<double>("pvErrorScaling");
00157 
00158         if (constraint == CONSTRAINT_PV_BEAMSPOT_SIZE ||
00159             constraint == CONSTRAINT_PV_BS_Z_ERRORS_SCALED ||
00160             constraint == CONSTRAINT_BEAMSPOT ||
00161             constraint == CONSTRAINT_PV_PRIMARIES_IN_FIT )
00162             beamSpotTag = params.getParameter<edm::InputTag>("beamSpotTag");
00163         useExternalSV = false;
00164         if(params.existsAs<bool>("useExternalSV")) useExternalSV = params.getParameter<bool> ("useExternalSV");
00165         if(useExternalSV) {
00166            extSVCollection  = params.getParameter<edm::InputTag>("extSVCollection");
00167            extSVDeltaRToJet = params.getParameter<double>("extSVDeltaRToJet");
00168         }
00169         produces<SecondaryVertexTagInfoCollection>();
00170 }
00171 
00172 SecondaryVertexProducer::~SecondaryVertexProducer()
00173 {
00174 }
00175 
00176 namespace {
00177         struct SVBuilder :
00178                 public std::unary_function<const Vertex&, SecondaryVertex> {
00179 
00180                 SVBuilder(const reco::Vertex &pv,
00181                           const GlobalVector &direction,
00182                           bool withPVError) :
00183                         pv(pv), direction(direction),
00184                         withPVError(withPVError) {}
00185 
00186                 SecondaryVertex operator () (const reco::Vertex &sv) const
00187                 { return SecondaryVertex(pv, sv, direction, withPVError); }
00188 
00189                 const Vertex            &pv;
00190                 const GlobalVector      &direction;
00191                 bool                    withPVError;
00192         };
00193 
00194         struct SVFilter :
00195                 public std::unary_function<const SecondaryVertex&, bool> {
00196 
00197                 SVFilter(const VertexFilter &filter, const Vertex &pv,
00198                          const GlobalVector &direction) :
00199                         filter(filter), pv(pv), direction(direction) {}
00200 
00201                 inline bool operator () (const SecondaryVertex &sv) const
00202                 { return !filter(pv, sv, direction); }
00203 
00204                 const VertexFilter      &filter;
00205                 const Vertex            &pv;
00206                 const GlobalVector      &direction;
00207         };
00208                         
00209 } // anonynmous namespace
00210 
00211 void SecondaryVertexProducer::produce(edm::Event &event,
00212                                       const edm::EventSetup &es)
00213 {
00214         typedef std::map<TrackBaseRef, TransientTrack,
00215                          RefToBaseLess<Track> > TransientTrackMap;
00216 
00217         edm::ESHandle<TransientTrackBuilder> trackBuilder;
00218         es.get<TransientTrackRecord>().get("TransientTrackBuilder",
00219                                            trackBuilder);
00220 
00221         edm::Handle<TrackIPTagInfoCollection> trackIPTagInfos;
00222         event.getByLabel(trackIPTagInfoLabel, trackIPTagInfos);
00223 
00224         // External Sec Vertex collection (e.g. for IVF usage)
00225         edm::Handle<reco::VertexCollection> extSecVertex;          
00226         if(useExternalSV) event.getByLabel(extSVCollection,extSecVertex);
00227                                                              
00228 
00229         edm::Handle<BeamSpot> beamSpot;
00230         unsigned int bsCovSrc[7] = { 0, };
00231         double sigmaZ = 0.0, beamWidth = 0.0;
00232         switch(constraint) {
00233             case CONSTRAINT_PV_BEAMSPOT_SIZE:
00234                 event.getByLabel(beamSpotTag,beamSpot);
00235                 bsCovSrc[3] = bsCovSrc[4] = bsCovSrc[5] = bsCovSrc[6] = 1;
00236                 sigmaZ = beamSpot->sigmaZ();
00237                 beamWidth = beamSpot->BeamWidthX();
00238                 break;
00239 
00240             case CONSTRAINT_PV_BS_Z_ERRORS_SCALED:
00241               event.getByLabel(beamSpotTag,beamSpot);
00242                 bsCovSrc[0] = bsCovSrc[1] = 2;
00243                 bsCovSrc[3] = bsCovSrc[4] = bsCovSrc[5] = 1;
00244                 sigmaZ = beamSpot->sigmaZ();
00245                 break;
00246 
00247             case CONSTRAINT_PV_ERROR_SCALED:
00248                 bsCovSrc[0] = bsCovSrc[1] = bsCovSrc[2] = 2;
00249                 break;
00250 
00251             case CONSTRAINT_BEAMSPOT:
00252             case CONSTRAINT_PV_PRIMARIES_IN_FIT:
00253                 event.getByLabel(beamSpotTag,beamSpot);
00254                 break;
00255 
00256             default:
00257                 /* nothing */;
00258         }
00259 
00260         std::auto_ptr<ConfigurableVertexReconstructor> vertexReco;
00261         std::auto_ptr<GhostTrackVertexFinder> vertexRecoGT;
00262         if (useGhostTrack)
00263                 vertexRecoGT.reset(new GhostTrackVertexFinder(
00264                         vtxRecoPSet.getParameter<double>("maxFitChi2"),
00265                         vtxRecoPSet.getParameter<double>("mergeThreshold"),
00266                         vtxRecoPSet.getParameter<double>("primcut"),
00267                         vtxRecoPSet.getParameter<double>("seccut"),
00268                         getGhostTrackFitType(vtxRecoPSet.getParameter<std::string>("fitType"))));
00269         else
00270                 vertexReco.reset(
00271                         new ConfigurableVertexReconstructor(vtxRecoPSet));
00272 
00273         TransientTrackMap primariesMap;
00274 
00275         // result secondary vertices
00276 
00277         std::auto_ptr<SecondaryVertexTagInfoCollection>
00278                         tagInfos(new SecondaryVertexTagInfoCollection);
00279 
00280         for(TrackIPTagInfoCollection::const_iterator iterJets =
00281                 trackIPTagInfos->begin(); iterJets != trackIPTagInfos->end();
00282                 ++iterJets) {
00283                 std::vector<SecondaryVertexTagInfo::IndexedTrackData> trackData;
00284 //                    std::cout << "Jet " << iterJets-trackIPTagInfos->begin() << std::endl; 
00285 
00286                 const Vertex &pv = *iterJets->primaryVertex();
00287 
00288                 std::set<TransientTrack> primaries;
00289                 if (constraint == CONSTRAINT_PV_PRIMARIES_IN_FIT) {
00290                         for(Vertex::trackRef_iterator iter = pv.tracks_begin();
00291                             iter != pv.tracks_end(); ++iter) {
00292                                 TransientTrackMap::iterator pos =
00293                                         primariesMap.lower_bound(*iter);
00294 
00295                                 if (pos != primariesMap.end() &&
00296                                     pos->first == *iter)
00297                                         primaries.insert(pos->second);
00298                                 else {
00299                                         TransientTrack track =
00300                                                 trackBuilder->build(
00301                                                         iter->castTo<TrackRef>());
00302                                         primariesMap.insert(pos,
00303                                                 std::make_pair(*iter, track));
00304                                         primaries.insert(track);
00305                                 }
00306                         }
00307                 }
00308 
00309                 edm::RefToBase<Jet> jetRef = iterJets->jet();
00310 
00311                 GlobalVector jetDir(jetRef->momentum().x(),
00312                                     jetRef->momentum().y(),
00313                                     jetRef->momentum().z());
00314 
00315                 std::vector<std::size_t> indices =
00316                                 iterJets->sortedIndexes(sortCriterium);
00317 
00318                 TrackRefVector trackRefs = iterJets->sortedTracks(indices);
00319 
00320                 const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
00321                                         iterJets->impactParameterData();
00322 
00323                 // build transient tracks used for vertex reconstruction
00324 
00325                 std::vector<TransientTrack> fitTracks;
00326                 std::vector<GhostTrackState> gtStates;
00327                 std::auto_ptr<GhostTrackPrediction> gtPred;
00328                 if (useGhostTrack)
00329                         gtPred.reset(new GhostTrackPrediction(
00330                                                 *iterJets->ghostTrack()));
00331 
00332                 for(unsigned int i = 0; i < indices.size(); i++) {
00333                         typedef SecondaryVertexTagInfo::IndexedTrackData IndexedTrackData;
00334 
00335                         const TrackRef &trackRef = trackRefs[i];
00336 
00337                         trackData.push_back(IndexedTrackData());
00338                         trackData.back().first = indices[i];
00339 
00340                         // select tracks for SV finder
00341 
00342                         if (!trackSelector(*trackRef, ipData[indices[i]], *jetRef,
00343                                            RecoVertex::convertPos(
00344                                                         pv.position()))) {
00345                                 trackData.back().second.svStatus =
00346                                         SecondaryVertexTagInfo::TrackData::trackSelected;
00347                                 continue;
00348                         }
00349 
00350                         TransientTrackMap::const_iterator pos =
00351                                         primariesMap.find(
00352                                                 TrackBaseRef(trackRef));
00353                         TransientTrack fitTrack;
00354                         if (pos != primariesMap.end()) {
00355                                 primaries.erase(pos->second);
00356                                 fitTrack = pos->second;
00357                         } else
00358                                 fitTrack = trackBuilder->build(trackRef);
00359                         fitTracks.push_back(fitTrack);
00360 
00361                         trackData.back().second.svStatus =
00362                                 SecondaryVertexTagInfo::TrackData::trackUsedForVertexFit;
00363 
00364                         if (useGhostTrack) {
00365                                 GhostTrackState gtState(fitTrack);
00366                                 GlobalPoint pos =
00367                                         ipData[indices[i]].closestToGhostTrack;
00368                                 gtState.linearize(*gtPred, true,
00369                                                   gtPred->lambda(pos));
00370                                 gtState.setWeight(ipData[indices[i]].ghostTrackWeight);
00371                                 gtStates.push_back(gtState);
00372                         }
00373                 }
00374 
00375                 std::auto_ptr<GhostTrack> ghostTrack;
00376                 if (useGhostTrack)
00377                         ghostTrack.reset(new GhostTrack(
00378                                 GhostTrackPrediction(
00379                                         RecoVertex::convertPos(pv.position()),
00380                                         RecoVertex::convertError(pv.error()),
00381                                         GlobalVector(
00382                                                 iterJets->ghostTrack()->px(),
00383                                                 iterJets->ghostTrack()->py(),
00384                                                 iterJets->ghostTrack()->pz()),
00385                                         0.05),
00386                                 *gtPred, gtStates,
00387                                 iterJets->ghostTrack()->chi2(),
00388                                 iterJets->ghostTrack()->ndof()));
00389 
00390                 // perform actual vertex finding
00391 
00392 
00393                 std::vector<reco::Vertex>       extAssoCollection;    
00394                 std::vector<TransientVertex> fittedSVs;
00395                 std::vector<SecondaryVertex> SVs;
00396                 if(!useExternalSV){ 
00397                   switch(constraint)   {
00398                     case CONSTRAINT_NONE:
00399                         if (useGhostTrack)
00400                                 fittedSVs = vertexRecoGT->vertices(
00401                                                 pv, *ghostTrack);
00402                         else
00403                                 fittedSVs = vertexReco->vertices(fitTracks);
00404                         break;
00405 
00406                     case CONSTRAINT_BEAMSPOT:
00407                         if (useGhostTrack)
00408                                 fittedSVs = vertexRecoGT->vertices(
00409                                                 pv, *beamSpot, *ghostTrack);
00410                         else
00411                                 fittedSVs = vertexReco->vertices(fitTracks,
00412                                                                  *beamSpot);
00413                         break;
00414 
00415                     case CONSTRAINT_PV_BEAMSPOT_SIZE:
00416                     case CONSTRAINT_PV_BS_Z_ERRORS_SCALED:
00417                     case CONSTRAINT_PV_ERROR_SCALED: {
00418                         BeamSpot::CovarianceMatrix cov;
00419                         for(unsigned int i = 0; i < 7; i++) {
00420                                 unsigned int covSrc = bsCovSrc[i];
00421                                 for(unsigned int j = 0; j < 7; j++) {
00422                                         double v=0.0;
00423                                         if (!covSrc || bsCovSrc[j] != covSrc)
00424                                                 v = 0.0;
00425                                         else if (covSrc == 1)
00426                                                 v = beamSpot->covariance(i, j);
00427                                         else if (j<3 && i<3)
00428                                                 v = pv.covariance(i, j) *
00429                                                     constraintScaling;
00430                                         cov(i, j) = v;
00431                                 }
00432                         }
00433 
00434                         BeamSpot bs(pv.position(), sigmaZ,
00435                                     beamSpot.isValid() ? beamSpot->dxdz() : 0.,
00436                                     beamSpot.isValid() ? beamSpot->dydz() : 0.,
00437                                     beamWidth, cov, BeamSpot::Unknown);
00438 
00439                         if (useGhostTrack)
00440                                 fittedSVs = vertexRecoGT->vertices(
00441                                                 pv, bs, *ghostTrack);
00442                         else
00443                                 fittedSVs = vertexReco->vertices(fitTracks, bs);
00444                     }   break;
00445 
00446                     case CONSTRAINT_PV_PRIMARIES_IN_FIT: {
00447                         std::vector<TransientTrack> primaries_(
00448                                         primaries.begin(), primaries.end());
00449                         if (useGhostTrack)
00450                                 fittedSVs = vertexRecoGT->vertices(
00451                                                 pv, *beamSpot, primaries_,
00452                                                 *ghostTrack);
00453                         else
00454                                 fittedSVs = vertexReco->vertices(
00455                                                 primaries_, fitTracks,
00456                                                 *beamSpot);
00457                     }   break;
00458                 }
00459 
00460                 // build combined SV information and filter
00461 
00462                 SVBuilder svBuilder(pv, jetDir, withPVError);
00463                 std::remove_copy_if(boost::make_transform_iterator(
00464                                         fittedSVs.begin(), svBuilder),
00465                                     boost::make_transform_iterator(
00466                                         fittedSVs.end(), svBuilder),
00467                                     std::back_inserter(SVs),
00468                                     SVFilter(vertexFilter, pv, jetDir));
00469 
00470                 // clean up now unneeded collections
00471              }else{
00472                   for(size_t iExtSv = 0; iExtSv < extSecVertex->size(); iExtSv++){
00473                       const reco::Vertex & extVertex = (*extSecVertex)[iExtSv];
00474 //                    GlobalVector vtxDir = GlobalVector(extVertex.p4().X(),extVertex.p4().Y(),extVertex.p4().Z());
00475 //                     if(Geom::deltaR(extVertex.position() - pv.position(), vtxDir)>0.2) continue; //pointing angle
00476 //                    std::cout << " dR " << iExtSv << " " << Geom::deltaR( ( extVertex.position() - pv.position() ), jetDir ) << "eta: " << ( extVertex.position() - pv.position()).eta() << " vs " << jetDir.eta() << " phi: "  << ( extVertex.position() - pv.position()).phi() << " vs  " << jetDir.phi() <<  std::endl; 
00477                       if( Geom::deltaR( ( extVertex.position() - pv.position() ), jetDir ) >  extSVDeltaRToJet || extVertex.p4().M() < 0.3)
00478                         continue;
00479 //                    std::cout << " SV added " << iExtSv << std::endl; 
00480                       extAssoCollection.push_back( extVertex );
00481                    }
00482                    SVBuilder svBuilder(pv, jetDir, withPVError);
00483                    std::remove_copy_if(boost::make_transform_iterator( extAssoCollection.begin(), svBuilder),
00484                                 boost::make_transform_iterator(extAssoCollection.end(), svBuilder),
00485                                 std::back_inserter(SVs),
00486                                 SVFilter(vertexFilter, pv, jetDir));
00487 
00488 
00489                  }
00490 //              std::cout << "size: " << SVs.size() << std::endl; 
00491                 gtPred.reset();
00492                 ghostTrack.reset();
00493                 gtStates.clear();
00494                 fitTracks.clear();
00495                 fittedSVs.clear();
00496                 extAssoCollection.clear();
00497 
00498                 // sort SVs by importance
00499 
00500                 std::vector<unsigned int> vtxIndices = vertexSorting(SVs);
00501 
00502                 std::vector<SecondaryVertexTagInfo::VertexData> svData;
00503 
00504                 svData.resize(vtxIndices.size());
00505                 for(unsigned int idx = 0; idx < vtxIndices.size(); idx++) {
00506                         const SecondaryVertex &sv = SVs[vtxIndices[idx]];
00507 
00508                         svData[idx].vertex = sv;
00509                         svData[idx].dist2d = sv.dist2d();
00510                         svData[idx].dist3d = sv.dist3d();
00511                         svData[idx].direction =
00512                                 GlobalVector(sv.x() - pv.x(),
00513                                              sv.y() - pv.y(),
00514                                              sv.z() - pv.z());
00515 
00516                         // mark tracks successfully used in vertex fit
00517 
00518                         for(Vertex::trackRef_iterator iter = sv.tracks_begin();
00519                             iter != sv.tracks_end(); iter++) {
00520                                 if (sv.trackWeight(*iter) < minTrackWeight)
00521                                         continue;
00522 
00523                                 TrackRefVector::const_iterator pos =
00524                                         std::find(trackRefs.begin(), trackRefs.end(),
00525                                                   iter->castTo<TrackRef>());
00526 
00527                                 if (pos == trackRefs.end() ) {
00528                                    if(!useExternalSV)
00529                                         throw cms::Exception("TrackNotFound")
00530                                                 << "Could not find track from secondary "
00531                                                    "vertex in original tracks."
00532                                                 << std::endl;
00533                                 } else {
00534                                 unsigned int index = pos - trackRefs.begin();
00535                                 trackData[index].second.svStatus =
00536                                         (SecondaryVertexTagInfo::TrackData::Status)
00537                                         ((unsigned int)SecondaryVertexTagInfo::TrackData::trackAssociatedToVertex + idx);
00538                                 }
00539                         }
00540                 }
00541 
00542                 // fill result into tag infos
00543 
00544                 tagInfos->push_back(
00545                         SecondaryVertexTagInfo(
00546                                 trackData, svData, SVs.size(),
00547                                 TrackIPTagInfoRef(trackIPTagInfos,
00548                                         iterJets - trackIPTagInfos->begin())));
00549         }
00550 
00551         event.put(tagInfos);
00552 }
00553 
00554 //define this as a plug-in
00555 DEFINE_FWK_MODULE(SecondaryVertexProducer);