![]() |
![]() |
Public Member Functions | |
virtual void | produce (edm::Event &event, const edm::EventSetup &es) |
SecondaryVertexProducer (const edm::ParameterSet ¶ms) | |
~SecondaryVertexProducer () | |
Private Attributes | |
double | minTrackWeight |
TrackIPTagInfo::SortCriteria | sortCriterium |
const edm::InputTag | trackIPTagInfoLabel |
TrackSelector | trackSelector |
bool | useBeamConstraint |
VertexFilter | vertexFilter |
VertexSorting | vertexSorting |
edm::ParameterSet | vtxRecoPSet |
bool | withPVError |
Definition at line 43 of file SecondaryVertexProducer.cc.
SecondaryVertexProducer::SecondaryVertexProducer | ( | const edm::ParameterSet & | params | ) | [explicit] |
Definition at line 62 of file SecondaryVertexProducer.cc.
00063 : 00064 trackIPTagInfoLabel(params.getParameter<edm::InputTag>("trackIPTagInfos")), 00065 sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))), 00066 trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")), 00067 useBeamConstraint(params.getParameter<bool>("useBeamConstraint")), 00068 vtxRecoPSet(params.getParameter<edm::ParameterSet>("vertexReco")), 00069 withPVError(params.getParameter<bool>("usePVError")), 00070 minTrackWeight(params.getParameter<double>("minimumTrackWeight")), 00071 vertexFilter(params.getParameter<edm::ParameterSet>("vertexCuts")), 00072 vertexSorting(params.getParameter<edm::ParameterSet>("vertexSelection")) 00073 { 00074 produces<SecondaryVertexTagInfoCollection>(); 00075 }
SecondaryVertexProducer::~SecondaryVertexProducer | ( | ) |
void SecondaryVertexProducer::produce | ( | edm::Event & | event, | |
const edm::EventSetup & | es | |||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 116 of file SecondaryVertexProducer.cc.
References edm::RefVector< C, T, F >::begin(), RecoVertex::convertPos(), reco::SecondaryVertex::dist2d(), reco::SecondaryVertex::dist3d(), edm::RefVector< C, T, F >::end(), lat::endl(), find(), edm::EventSetup::get(), i, index, iter, jetRef, minTrackWeight, reco::Vertex::position(), edm::RefVector< C, T, F >::push_back(), pv, sortCriterium, sv, tagInfos, trackIPTagInfoLabel, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), trackSelector, reco::Vertex::trackWeight(), useBeamConstraint, vertexFilter, vertexSorting, ConfigurableVertexReconstructor::vertices(), vtxRecoPSet, withPVError, reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().
00118 { 00119 edm::ESHandle<TransientTrackBuilder> trackBuilder; 00120 es.get<TransientTrackRecord>().get("TransientTrackBuilder", 00121 trackBuilder); 00122 00123 edm::Handle<TrackIPTagInfoCollection> trackIPTagInfos; 00124 event.getByLabel(trackIPTagInfoLabel, trackIPTagInfos); 00125 00126 edm::Handle<BeamSpot> beamSpot; 00127 if (useBeamConstraint) 00128 event.getByType(beamSpot); 00129 00130 ConfigurableVertexReconstructor vertexReco(vtxRecoPSet); 00131 00132 // result secondary vertices 00133 00134 std::auto_ptr<SecondaryVertexTagInfoCollection> 00135 tagInfos(new SecondaryVertexTagInfoCollection); 00136 00137 for(TrackIPTagInfoCollection::const_iterator iterJets = 00138 trackIPTagInfos->begin(); iterJets != trackIPTagInfos->end(); 00139 ++iterJets) { 00140 00141 std::vector<SecondaryVertexTagInfo::IndexedTrackData> trackData; 00142 00143 const Vertex &pv = *iterJets->primaryVertex(); 00144 00145 edm::RefToBase<Jet> jetRef = iterJets->jet(); 00146 00147 GlobalVector jetDir(jetRef->momentum().x(), 00148 jetRef->momentum().y(), 00149 jetRef->momentum().z()); 00150 00151 std::vector<std::size_t> indices = 00152 iterJets->sortedIndexes(sortCriterium); 00153 00154 TrackRefVector trackRefs = iterJets->sortedTracks(indices); 00155 00156 const std::vector<TrackIPTagInfo::TrackIPData> &ipData = 00157 iterJets->impactParameterData(); 00158 00159 // build transient tracks used for vertex reconstruction 00160 00161 std::vector<TransientTrack> fitTracks; 00162 for(unsigned int i = 0; i < indices.size(); i++) { 00163 typedef SecondaryVertexTagInfo::IndexedTrackData IndexedTrackData; 00164 00165 const TrackRef &trackRef = trackRefs[i]; 00166 00167 trackData.push_back(IndexedTrackData()); 00168 trackData.back().first = indices[i]; 00169 00170 // select tracks for SV finder 00171 00172 if (trackSelector(*trackRef, ipData[i], *jetRef, 00173 RecoVertex::convertPos( 00174 pv.position()))) { 00175 fitTracks.push_back( 00176 trackBuilder->build(trackRef)); 00177 trackData.back().second.svStatus = 00178 SecondaryVertexTagInfo::TrackData::trackUsedForVertexFit; 00179 } else 00180 trackData.back().second.svStatus = 00181 SecondaryVertexTagInfo::TrackData::trackSelected; 00182 } 00183 00184 // perform actual vertex finding 00185 00186 std::vector<TransientVertex> fittedSVs; 00187 if (useBeamConstraint) 00188 fittedSVs = vertexReco.vertices(fitTracks, *beamSpot); 00189 else 00190 fittedSVs = vertexReco.vertices(fitTracks); 00191 00192 // build combined SV information and filter 00193 00194 std::vector<SecondaryVertex> SVs; 00195 SVBuilder svBuilder(pv, jetDir, withPVError); 00196 std::remove_copy_if(boost::make_transform_iterator( 00197 fittedSVs.begin(), svBuilder), 00198 boost::make_transform_iterator( 00199 fittedSVs.end(), svBuilder), 00200 std::back_inserter(SVs), 00201 SVFilter(vertexFilter, pv, jetDir)); 00202 00203 // clean up now unneeded collections 00204 00205 fitTracks.clear(); 00206 fittedSVs.clear(); 00207 00208 // sort SVs by importance 00209 00210 std::vector<unsigned int> vtxIndices = vertexSorting(SVs); 00211 00212 std::vector<SecondaryVertexTagInfo::VertexData> svData; 00213 00214 svData.resize(vtxIndices.size()); 00215 for(unsigned int idx = 0; idx < vtxIndices.size(); idx++) { 00216 const SecondaryVertex &sv = SVs[vtxIndices[idx]]; 00217 00218 svData[idx].vertex = sv; 00219 svData[idx].dist2d = sv.dist2d(); 00220 svData[idx].dist3d = sv.dist3d(); 00221 svData[idx].direction = 00222 GlobalVector(sv.x() - pv.x(), 00223 sv.y() - pv.y(), 00224 sv.z() - pv.z()); 00225 00226 // mark tracks successfully used in vertex fit 00227 00228 for(Vertex::trackRef_iterator iter = sv.tracks_begin(); 00229 iter != sv.tracks_end(); iter++) { 00230 if (sv.trackWeight(*iter) < minTrackWeight) 00231 continue; 00232 00233 TrackRefVector::const_iterator pos = 00234 std::find(trackRefs.begin(), trackRefs.end(), 00235 iter->castTo<TrackRef>()); 00236 if (pos == trackRefs.end()) 00237 throw cms::Exception("TrackNotFound") 00238 << "Could not find track from secondary " 00239 "vertex in original tracks." 00240 << std::endl; 00241 00242 unsigned int index = pos - trackRefs.begin(); 00243 trackData[index].second.svStatus = 00244 (SecondaryVertexTagInfo::TrackData::Status) 00245 ((unsigned int)SecondaryVertexTagInfo::TrackData::trackAssociatedToVertex + idx); 00246 } 00247 } 00248 00249 // fill result into tag infos 00250 00251 tagInfos->push_back( 00252 SecondaryVertexTagInfo( 00253 trackData, svData, SVs.size(), 00254 TrackIPTagInfoRef(trackIPTagInfos, 00255 iterJets - trackIPTagInfos->begin()))); 00256 } 00257 00258 event.put(tagInfos); 00259 }
double SecondaryVertexProducer::minTrackWeight [private] |
const edm::InputTag SecondaryVertexProducer::trackIPTagInfoLabel [private] |
bool SecondaryVertexProducer::withPVError [private] |