Public Member Functions | |
virtual void | produce (edm::Event &event, const edm::EventSetup &es) |
SecondaryVertexProducer (const edm::ParameterSet ¶ms) | |
~SecondaryVertexProducer () | |
Private Types | |
enum | ConstraintType { CONSTRAINT_NONE = 0, CONSTRAINT_BEAMSPOT, CONSTRAINT_PV_BEAMSPOT_SIZE, CONSTRAINT_PV_BS_Z_ERRORS_SCALED, CONSTRAINT_PV_ERROR_SCALED, CONSTRAINT_PV_PRIMARIES_IN_FIT } |
Static Private Member Functions | |
static ConstraintType | getConstraintType (const std::string &name) |
Private Attributes | |
edm::InputTag | beamSpotTag |
ConstraintType | constraint |
double | constraintScaling |
double | minTrackWeight |
TrackIPTagInfo::SortCriteria | sortCriterium |
const edm::InputTag | trackIPTagInfoLabel |
TrackSelector | trackSelector |
bool | useGhostTrack |
VertexFilter | vertexFilter |
VertexSorting | vertexSorting |
edm::ParameterSet | vtxRecoPSet |
bool | withPVError |
Definition at line 63 of file SecondaryVertexProducer.cc.
enum SecondaryVertexProducer::ConstraintType [private] |
CONSTRAINT_NONE | |
CONSTRAINT_BEAMSPOT | |
CONSTRAINT_PV_BEAMSPOT_SIZE | |
CONSTRAINT_PV_BS_Z_ERRORS_SCALED | |
CONSTRAINT_PV_ERROR_SCALED | |
CONSTRAINT_PV_PRIMARIES_IN_FIT |
Definition at line 71 of file SecondaryVertexProducer.cc.
SecondaryVertexProducer::SecondaryVertexProducer | ( | const edm::ParameterSet & | params | ) | [explicit] |
Definition at line 135 of file SecondaryVertexProducer.cc.
References beamSpotTag, constraint, CONSTRAINT_BEAMSPOT, CONSTRAINT_PV_BEAMSPOT_SIZE, CONSTRAINT_PV_BS_Z_ERRORS_SCALED, CONSTRAINT_PV_ERROR_SCALED, CONSTRAINT_PV_PRIMARIES_IN_FIT, constraintScaling, and edm::ParameterSet::getParameter().
: trackIPTagInfoLabel(params.getParameter<edm::InputTag>("trackIPTagInfos")), sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))), trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")), constraint(getConstraintType(params.getParameter<std::string>("constraint"))), constraintScaling(1.0), vtxRecoPSet(params.getParameter<edm::ParameterSet>("vertexReco")), useGhostTrack(vtxRecoPSet.getParameter<std::string>("finder") == "gtvr"), withPVError(params.getParameter<bool>("usePVError")), minTrackWeight(params.getParameter<double>("minimumTrackWeight")), vertexFilter(params.getParameter<edm::ParameterSet>("vertexCuts")), vertexSorting(params.getParameter<edm::ParameterSet>("vertexSelection")) { if (constraint == CONSTRAINT_PV_ERROR_SCALED || constraint == CONSTRAINT_PV_BS_Z_ERRORS_SCALED) constraintScaling = params.getParameter<double>("pvErrorScaling"); if (constraint == CONSTRAINT_PV_BEAMSPOT_SIZE || constraint == CONSTRAINT_PV_BS_Z_ERRORS_SCALED || constraint == CONSTRAINT_BEAMSPOT || constraint == CONSTRAINT_PV_PRIMARIES_IN_FIT ) beamSpotTag = params.getParameter<edm::InputTag>("beamSpotTag"); produces<SecondaryVertexTagInfoCollection>(); }
SecondaryVertexProducer::~SecondaryVertexProducer | ( | ) |
Definition at line 162 of file SecondaryVertexProducer.cc.
{ }
SecondaryVertexProducer::ConstraintType SecondaryVertexProducer::getConstraintType | ( | const std::string & | name | ) | [static, private] |
Definition at line 97 of file SecondaryVertexProducer.cc.
References Exception.
{ if (name == "None") return CONSTRAINT_NONE; else if (name == "BeamSpot") return CONSTRAINT_BEAMSPOT; else if (name == "BeamSpot+PVPosition") return CONSTRAINT_PV_BEAMSPOT_SIZE; else if (name == "BeamSpotZ+PVErrorScaledXY") return CONSTRAINT_PV_BS_Z_ERRORS_SCALED; else if (name == "PVErrorScaled") return CONSTRAINT_PV_ERROR_SCALED; else if (name == "BeamSpot+PVTracksInFit") return CONSTRAINT_PV_PRIMARIES_IN_FIT; else throw cms::Exception("InvalidArgument") << "SecondaryVertexProducer: ``constraint'' parameter " "value \"" << name << "\" not understood." << std::endl; }
void SecondaryVertexProducer::produce | ( | edm::Event & | event, |
const edm::EventSetup & | es | ||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 201 of file SecondaryVertexProducer.cc.
References beamSpotTag, edm::RefVector< C, T, F >::begin(), constraint, CONSTRAINT_BEAMSPOT, CONSTRAINT_NONE, CONSTRAINT_PV_BEAMSPOT_SIZE, CONSTRAINT_PV_BS_Z_ERRORS_SCALED, CONSTRAINT_PV_ERROR_SCALED, CONSTRAINT_PV_PRIMARIES_IN_FIT, constraintScaling, RecoVertex::convertError(), RecoVertex::convertPos(), reco::Vertex::covariance(), reco::SecondaryVertex::dist2d(), reco::SecondaryVertex::dist3d(), edm::RefVector< C, T, F >::end(), reco::Vertex::error(), spr::find(), edm::EventSetup::get(), getGhostTrackFitType(), edm::ParameterSet::getParameter(), ghostTrackES_cfi::ghostTrack, i, getHLTprescales::index, edm::HandleBase::isValid(), j, reco::GhostTrackState::linearize(), python::multivaluedict::map(), minTrackWeight, pos, reco::Vertex::position(), edm::RefVector< C, T, F >::push_back(), reco::GhostTrackState::setWeight(), sortCriterium, step1_ZMM_7Tev::Status, trackIPTagInfoLabel, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), trackSelector, reco::Vertex::trackWeight(), Unknown, useGhostTrack, v, vertexFilter, vertexSorting, vtxRecoPSet, withPVError, reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().
{ typedef std::map<TrackBaseRef, TransientTrack, RefToBaseLess<Track> > TransientTrackMap; edm::ESHandle<TransientTrackBuilder> trackBuilder; es.get<TransientTrackRecord>().get("TransientTrackBuilder", trackBuilder); edm::Handle<TrackIPTagInfoCollection> trackIPTagInfos; event.getByLabel(trackIPTagInfoLabel, trackIPTagInfos); edm::Handle<BeamSpot> beamSpot; unsigned int bsCovSrc[7] = { 0, }; double sigmaZ = 0.0, beamWidth = 0.0; switch(constraint) { case CONSTRAINT_PV_BEAMSPOT_SIZE: event.getByLabel(beamSpotTag,beamSpot); bsCovSrc[3] = bsCovSrc[4] = bsCovSrc[5] = bsCovSrc[6] = 1; sigmaZ = beamSpot->sigmaZ(); beamWidth = beamSpot->BeamWidthX(); break; case CONSTRAINT_PV_BS_Z_ERRORS_SCALED: event.getByLabel(beamSpotTag,beamSpot); bsCovSrc[0] = bsCovSrc[1] = 2; bsCovSrc[3] = bsCovSrc[4] = bsCovSrc[5] = 1; sigmaZ = beamSpot->sigmaZ(); break; case CONSTRAINT_PV_ERROR_SCALED: bsCovSrc[0] = bsCovSrc[1] = bsCovSrc[2] = 2; break; case CONSTRAINT_BEAMSPOT: case CONSTRAINT_PV_PRIMARIES_IN_FIT: event.getByLabel(beamSpotTag,beamSpot); break; default: /* nothing */; } std::auto_ptr<ConfigurableVertexReconstructor> vertexReco; std::auto_ptr<GhostTrackVertexFinder> vertexRecoGT; if (useGhostTrack) vertexRecoGT.reset(new GhostTrackVertexFinder( vtxRecoPSet.getParameter<double>("maxFitChi2"), vtxRecoPSet.getParameter<double>("mergeThreshold"), vtxRecoPSet.getParameter<double>("primcut"), vtxRecoPSet.getParameter<double>("seccut"), getGhostTrackFitType(vtxRecoPSet.getParameter<std::string>("fitType")))); else vertexReco.reset( new ConfigurableVertexReconstructor(vtxRecoPSet)); TransientTrackMap primariesMap; // result secondary vertices std::auto_ptr<SecondaryVertexTagInfoCollection> tagInfos(new SecondaryVertexTagInfoCollection); for(TrackIPTagInfoCollection::const_iterator iterJets = trackIPTagInfos->begin(); iterJets != trackIPTagInfos->end(); ++iterJets) { std::vector<SecondaryVertexTagInfo::IndexedTrackData> trackData; const Vertex &pv = *iterJets->primaryVertex(); std::set<TransientTrack> primaries; if (constraint == CONSTRAINT_PV_PRIMARIES_IN_FIT) { for(Vertex::trackRef_iterator iter = pv.tracks_begin(); iter != pv.tracks_end(); ++iter) { TransientTrackMap::iterator pos = primariesMap.lower_bound(*iter); if (pos != primariesMap.end() && pos->first == *iter) primaries.insert(pos->second); else { TransientTrack track = trackBuilder->build( iter->castTo<TrackRef>()); primariesMap.insert(pos, std::make_pair(*iter, track)); primaries.insert(track); } } } edm::RefToBase<Jet> jetRef = iterJets->jet(); GlobalVector jetDir(jetRef->momentum().x(), jetRef->momentum().y(), jetRef->momentum().z()); std::vector<std::size_t> indices = iterJets->sortedIndexes(sortCriterium); TrackRefVector trackRefs = iterJets->sortedTracks(indices); const std::vector<TrackIPTagInfo::TrackIPData> &ipData = iterJets->impactParameterData(); // build transient tracks used for vertex reconstruction std::vector<TransientTrack> fitTracks; std::vector<GhostTrackState> gtStates; std::auto_ptr<GhostTrackPrediction> gtPred; if (useGhostTrack) gtPred.reset(new GhostTrackPrediction( *iterJets->ghostTrack())); for(unsigned int i = 0; i < indices.size(); i++) { typedef SecondaryVertexTagInfo::IndexedTrackData IndexedTrackData; const TrackRef &trackRef = trackRefs[i]; trackData.push_back(IndexedTrackData()); trackData.back().first = indices[i]; // select tracks for SV finder if (!trackSelector(*trackRef, ipData[indices[i]], *jetRef, RecoVertex::convertPos( pv.position()))) { trackData.back().second.svStatus = SecondaryVertexTagInfo::TrackData::trackSelected; continue; } TransientTrackMap::const_iterator pos = primariesMap.find( TrackBaseRef(trackRef)); TransientTrack fitTrack; if (pos != primariesMap.end()) { primaries.erase(pos->second); fitTrack = pos->second; } else fitTrack = trackBuilder->build(trackRef); fitTracks.push_back(fitTrack); trackData.back().second.svStatus = SecondaryVertexTagInfo::TrackData::trackUsedForVertexFit; if (useGhostTrack) { GhostTrackState gtState(fitTrack); GlobalPoint pos = ipData[indices[i]].closestToGhostTrack; gtState.linearize(*gtPred, true, gtPred->lambda(pos)); gtState.setWeight(ipData[indices[i]].ghostTrackWeight); gtStates.push_back(gtState); } } std::auto_ptr<GhostTrack> ghostTrack; if (useGhostTrack) ghostTrack.reset(new GhostTrack( GhostTrackPrediction( RecoVertex::convertPos(pv.position()), RecoVertex::convertError(pv.error()), GlobalVector( iterJets->ghostTrack()->px(), iterJets->ghostTrack()->py(), iterJets->ghostTrack()->pz()), 0.05), *gtPred, gtStates, iterJets->ghostTrack()->chi2(), iterJets->ghostTrack()->ndof())); // perform actual vertex finding std::vector<TransientVertex> fittedSVs; switch(constraint) { case CONSTRAINT_NONE: if (useGhostTrack) fittedSVs = vertexRecoGT->vertices( pv, *ghostTrack); else fittedSVs = vertexReco->vertices(fitTracks); break; case CONSTRAINT_BEAMSPOT: if (useGhostTrack) fittedSVs = vertexRecoGT->vertices( pv, *beamSpot, *ghostTrack); else fittedSVs = vertexReco->vertices(fitTracks, *beamSpot); break; case CONSTRAINT_PV_BEAMSPOT_SIZE: case CONSTRAINT_PV_BS_Z_ERRORS_SCALED: case CONSTRAINT_PV_ERROR_SCALED: { BeamSpot::CovarianceMatrix cov; for(unsigned int i = 0; i < 7; i++) { unsigned int covSrc = bsCovSrc[i]; for(unsigned int j = 0; j < 7; j++) { double v=0.0; if (!covSrc || bsCovSrc[j] != covSrc) v = 0.0; else if (covSrc == 1) v = beamSpot->covariance(i, j); else if (j<3 && i<3) v = pv.covariance(i, j) * constraintScaling; cov(i, j) = v; } } BeamSpot bs(pv.position(), sigmaZ, beamSpot.isValid() ? beamSpot->dxdz() : 0., beamSpot.isValid() ? beamSpot->dydz() : 0., beamWidth, cov, BeamSpot::Unknown); if (useGhostTrack) fittedSVs = vertexRecoGT->vertices( pv, bs, *ghostTrack); else fittedSVs = vertexReco->vertices(fitTracks, bs); } break; case CONSTRAINT_PV_PRIMARIES_IN_FIT: { std::vector<TransientTrack> primaries_( primaries.begin(), primaries.end()); if (useGhostTrack) fittedSVs = vertexRecoGT->vertices( pv, *beamSpot, primaries_, *ghostTrack); else fittedSVs = vertexReco->vertices( primaries_, fitTracks, *beamSpot); } break; } // build combined SV information and filter std::vector<SecondaryVertex> SVs; SVBuilder svBuilder(pv, jetDir, withPVError); std::remove_copy_if(boost::make_transform_iterator( fittedSVs.begin(), svBuilder), boost::make_transform_iterator( fittedSVs.end(), svBuilder), std::back_inserter(SVs), SVFilter(vertexFilter, pv, jetDir)); // clean up now unneeded collections gtPred.reset(); ghostTrack.reset(); gtStates.clear(); fitTracks.clear(); fittedSVs.clear(); // sort SVs by importance std::vector<unsigned int> vtxIndices = vertexSorting(SVs); std::vector<SecondaryVertexTagInfo::VertexData> svData; svData.resize(vtxIndices.size()); for(unsigned int idx = 0; idx < vtxIndices.size(); idx++) { const SecondaryVertex &sv = SVs[vtxIndices[idx]]; svData[idx].vertex = sv; svData[idx].dist2d = sv.dist2d(); svData[idx].dist3d = sv.dist3d(); svData[idx].direction = GlobalVector(sv.x() - pv.x(), sv.y() - pv.y(), sv.z() - pv.z()); // mark tracks successfully used in vertex fit for(Vertex::trackRef_iterator iter = sv.tracks_begin(); iter != sv.tracks_end(); iter++) { if (sv.trackWeight(*iter) < minTrackWeight) continue; TrackRefVector::const_iterator pos = std::find(trackRefs.begin(), trackRefs.end(), iter->castTo<TrackRef>()); if (pos == trackRefs.end()) throw cms::Exception("TrackNotFound") << "Could not find track from secondary " "vertex in original tracks." << std::endl; unsigned int index = pos - trackRefs.begin(); trackData[index].second.svStatus = (SecondaryVertexTagInfo::TrackData::Status) ((unsigned int)SecondaryVertexTagInfo::TrackData::trackAssociatedToVertex + idx); } } // fill result into tag infos tagInfos->push_back( SecondaryVertexTagInfo( trackData, svData, SVs.size(), TrackIPTagInfoRef(trackIPTagInfos, iterJets - trackIPTagInfos->begin()))); } event.put(tagInfos); }
Definition at line 83 of file SecondaryVertexProducer.cc.
Referenced by produce(), and SecondaryVertexProducer().
Definition at line 86 of file SecondaryVertexProducer.cc.
Referenced by produce(), and SecondaryVertexProducer().
double SecondaryVertexProducer::constraintScaling [private] |
Definition at line 87 of file SecondaryVertexProducer.cc.
Referenced by produce(), and SecondaryVertexProducer().
double SecondaryVertexProducer::minTrackWeight [private] |
Definition at line 91 of file SecondaryVertexProducer.cc.
Referenced by produce().
Definition at line 84 of file SecondaryVertexProducer.cc.
Referenced by produce().
const edm::InputTag SecondaryVertexProducer::trackIPTagInfoLabel [private] |
Definition at line 82 of file SecondaryVertexProducer.cc.
Referenced by produce().
Definition at line 85 of file SecondaryVertexProducer.cc.
Referenced by produce().
bool SecondaryVertexProducer::useGhostTrack [private] |
Definition at line 89 of file SecondaryVertexProducer.cc.
Referenced by produce().
Definition at line 92 of file SecondaryVertexProducer.cc.
Referenced by produce().
Definition at line 93 of file SecondaryVertexProducer.cc.
Referenced by produce().
Definition at line 88 of file SecondaryVertexProducer.cc.
Referenced by produce().
bool SecondaryVertexProducer::withPVError [private] |
Definition at line 90 of file SecondaryVertexProducer.cc.
Referenced by produce().