CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
TrackVertexArbitration< VTX > Class Template Reference

#include <TrackVertexArbitration.h>

Public Member Functions

 TrackVertexArbitration (const edm::ParameterSet &params)
 
std::vector< VTX > trackVertexArbitrator (edm::Handle< reco::BeamSpot > &beamSpot, const reco::Vertex &pv, std::vector< reco::TransientTrack > &selectedTracks, std::vector< VTX > &secondaryVertices)
 

Private Member Functions

bool trackFilterArbitrator (const reco::TransientTrack &track) const
 

Private Attributes

edm::InputTag beamSpotCollection
 
double distCut
 
double dLenFraction
 
double dRCut
 
double fitterRatio
 
double fitterSigmacut
 
double fitterTini
 
double maxTimeSignificance
 
edm::InputTag primaryVertexCollection
 
edm::InputTag secondaryVertexCollection
 
double sigCut
 
double sign
 
edm::InputTag trackCollection
 
int trackMinLayers
 
int trackMinPixels
 
double trackMinPt
 

Detailed Description

template<class VTX>
class TrackVertexArbitration< VTX >

Definition at line 58 of file TrackVertexArbitration.h.

Constructor & Destructor Documentation

template<class VTX >
TrackVertexArbitration< VTX >::TrackVertexArbitration ( const edm::ParameterSet params)

Definition at line 94 of file TrackVertexArbitration.h.

References TrackVertexArbitration< VTX >::dRCut, and TrackVertexArbitration< VTX >::sign.

94  :
95  primaryVertexCollection (params.getParameter<edm::InputTag>("primaryVertices")),
96  secondaryVertexCollection (params.getParameter<edm::InputTag>("secondaryVertices")),
97  trackCollection (params.getParameter<edm::InputTag>("tracks")),
98  beamSpotCollection (params.getParameter<edm::InputTag>("beamSpot")),
99  dRCut (params.getParameter<double>("dRCut")),
100  distCut (params.getParameter<double>("distCut")),
101  sigCut (params.getParameter<double>("sigCut")),
102  dLenFraction (params.getParameter<double>("dLenFraction")),
103  fitterSigmacut (params.getParameter<double>("fitterSigmacut")),
104  fitterTini (params.getParameter<double>("fitterTini")),
105  fitterRatio (params.getParameter<double>("fitterRatio")),
106  trackMinLayers (params.getParameter<int32_t>("trackMinLayers")),
107  trackMinPt (params.getParameter<double>("trackMinPt")),
108  trackMinPixels (params.getParameter<int32_t>("trackMinPixels")),
109  maxTimeSignificance (params.getParameter<double>("maxTimeSignificance"))
110 {
111  sign = 1.0;
112  if(dRCut < 0){sign = -1.0;}
113  dRCut*=dRCut;
114 }
T getParameter(std::string const &) const
edm::InputTag secondaryVertexCollection

Member Function Documentation

template<class VTX >
bool TrackVertexArbitration< VTX >::trackFilterArbitrator ( const reco::TransientTrack track) const
private

Definition at line 116 of file TrackVertexArbitration.h.

References reco::TrackBase::hitPattern(), reco::TransientTrack::isValid(), reco::HitPattern::numberOfValidPixelHits(), reco::TrackBase::pt(), reco::TransientTrack::track(), reco::HitPattern::trackerLayersWithMeasurement(), TrackVertexArbitration< VTX >::trackMinLayers, TrackVertexArbitration< VTX >::trackMinPixels, and TrackVertexArbitration< VTX >::trackMinPt.

Referenced by TrackVertexArbitration< VTX >::trackVertexArbitrator().

117 {
118  if(!track.isValid())
119  return false;
121  return false;
122  if (track.track().pt() < trackMinPt)
123  return false;
125  return false;
126 
127  return true;
128 }
bool isValid() const
Make the ReferenceCountingProxy method to check validity public.
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:557
double pt() const
track transverse momentum
Definition: TrackBase.h:654
const Track & track() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:479
int numberOfValidPixelHits() const
Definition: HitPattern.h:908
template<class VTX >
std::vector< VTX > TrackVertexArbitration< VTX >::trackVertexArbitrator ( edm::Handle< reco::BeamSpot > &  beamSpot,
const reco::Vertex pv,
std::vector< reco::TransientTrack > &  selectedTracks,
std::vector< VTX > &  secondaryVertices 
)

Definition at line 153 of file TrackVertexArbitration.h.

References funct::abs(), IPTools::absoluteImpactParameter3D(), RecoVertex::convertError(), RecoVertex::convertPos(), gather_cfg::cout, svhelper::cov33(), reco::Vertex::covariance(), deltaR(), reco::deltaR2(), VertexDistance3D::distance(), TrackVertexArbitration< VTX >::distCut, TrackVertexArbitration< VTX >::dLenFraction, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, TrackVertexArbitration< VTX >::dRCut, reco::TransientTrack::dtErrorExt(), reco::TransientTrack::field(), TrackVertexArbitration< VTX >::fitterRatio, TrackVertexArbitration< VTX >::fitterSigmacut, TrackVertexArbitration< VTX >::fitterTini, reco::TransientTrack::impactPointState(), edm::isFinite(), TransientVertex::isValid(), TrackVertexArbitration< VTX >::maxTimeSignificance, reco::Vertex::position(), funct::pow(), trackingPlots::selectedTracks, reco::TransientTrack::setBeamSpot(), TrackVertexArbitration< VTX >::sigCut, TrackVertexArbitration< VTX >::sign, Measurement1D::significance(), mathSSE::sqrt(), pfDeepBoostedJetPreprocessParams_cfi::sv, reco::TransientTrack::timeExt(), reco::TransientTrack::track(), TrackVertexArbitration< VTX >::trackFilterArbitrator(), trackWeight(), groupFilesInBlocks::tt, svhelper::updateVertexTime(), Measurement1D::value(), AdaptiveVertexFitter::vertex(), and w.

158 {
159  using namespace reco;
160 
161  //std::cout << "PV: " << pv.position() << std::endl;
162  VertexDistance3D dist;
163  AdaptiveVertexFitter theAdaptiveFitter(
169 
170 
171 
172  std::vector<VTX> recoVertices;
173  VertexDistance3D vdist;
174  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
175 
176  std::unordered_map<unsigned int, Measurement1D> cachedIP;
177  for(typename std::vector<VTX>::const_iterator sv = secondaryVertices.begin();
178  sv != secondaryVertices.end(); ++sv) {
179 
180  const double svTime(sv->t()), svTimeCov(svhelper::cov33(*sv));
181  const bool svHasTime = svTimeCov > 0.;
182 
183  GlobalPoint ssv(sv->position().x(),sv->position().y(),sv->position().z());
184  GlobalVector flightDir = ssv-ppv;
185 // std::cout << "Vertex : " << sv-secondaryVertices->begin() << " " << sv->position() << std::endl;
187  std::vector<reco::TransientTrack> selTracks;
188  for(unsigned int itrack = 0; itrack < selectedTracks.size(); itrack++){
189  TransientTrack & tt = (selectedTracks)[itrack];
190  if (!trackFilterArbitrator(tt)) continue;
191  tt.setBeamSpot(*beamSpot);
192  float w = trackWeight(*sv,tt);
193  Measurement1D ipv;
194  if( cachedIP.count(itrack) ) {
195  ipv=cachedIP[itrack];
196  } else {
197  std::pair<bool,Measurement1D> ipvp = IPTools::absoluteImpactParameter3D(tt,pv);
198  cachedIP[itrack]=ipvp.second;
199  ipv=ipvp.second;
200  }
201 
202  AnalyticalImpactPointExtrapolator extrapolator(tt.field());
203  TrajectoryStateOnSurface tsos = extrapolator.extrapolate(tt.impactPointState(), RecoVertex::convertPos(sv->position()));
204  if(! tsos.isValid()) continue;
205  GlobalPoint refPoint = tsos.globalPosition();
206  GlobalError refPointErr = tsos.cartesianError().position();
207  Measurement1D isv = dist.distance(VertexState(RecoVertex::convertPos(sv->position()),RecoVertex::convertError(sv->error())),VertexState(refPoint, refPointErr));
208 
209  float dR = Geom::deltaR2(( (sign > 0) ? flightDir : -flightDir),tt.track()); //.eta(), flightDir.phi(), tt.track().eta(), tt.track().phi());
210 
211  double timeSig = 0.;
212  if( svHasTime && edm::isFinite(tt.timeExt()) ) {
213  double tError = std::sqrt( std::pow( tt.dtErrorExt(), 2 ) + svTimeCov );
214  timeSig = std::abs(tt.timeExt() - svTime)/tError;
215  }
216 
217  if( w > 0 || ( isv.significance() < sigCut && isv.value() < distCut && isv.value() < dlen.value()*dLenFraction && timeSig < maxTimeSignificance) )
218  {
219 
220  if(( isv.value() < ipv.value() ) && isv.value() < distCut && isv.value() < dlen.value()*dLenFraction
221  && dR < dRCut && timeSig < maxTimeSignificance )
222  {
223 #ifdef VTXDEBUG
224  if(w > 0.5) std::cout << " = ";
225  else std::cout << " + ";
226 #endif
227  selTracks.push_back(tt);
228  } else
229  {
230 #ifdef VTXDEBUG
231  if(w > 0.5 && isv.value() > ipv.value() ) std::cout << " - ";
232  else std::cout << " ";
233 #endif
234  //add also the tracks used in previous fitting that are still closer to Sv than Pv
235  if(w > 0.5 && isv.value() <= ipv.value() && dR < dRCut && timeSig < maxTimeSignificance) {
236  selTracks.push_back(tt);
237 #ifdef VTXDEBUG
238  std::cout << " = ";
239 #endif
240  }
241  if(w > 0.5 && isv.value() <= ipv.value() && dR >= dRCut) {
242 #ifdef VTXDEBUG
243  std::cout << " - ";
244 #endif
245 
246  }
247 
248 
249  }
250 #ifdef VTXDEBUG
251 
252  std::cout << "t : " << itrack << " ref " << ref.key() << " w: " << w
253  << " svip: " << isv.significance() << " " << isv.value()
254  << " pvip: " << ipv.significance() << " " << ipv.value() << " res " << tt.track().residualX(0) << "," << tt.track().residualY(0)
255 // << " tpvip: " << itpv.second.significance() << " " << itpv.second.value() << " dr: " << dR << std::endl;
256 #endif
257 
258  }
259  else
260  {
261 #ifdef VTXDEBUG
262 
263  std::cout << " . t : " << itrack << " ref " << ref.key()<< " w: " << w
264  << " svip: " << isv.second.significance() << " " << isv.second.value()
265  << " pvip: " << ipv.significance() << " " << ipv.value() << " dr: " << dR << std::endl;
266 #endif
267 
268  }
269  }
270 
271  if(selTracks.size() >= 2)
272  {
273  TransientVertex singleFitVertex;
274  singleFitVertex = theAdaptiveFitter.vertex(selTracks,ssv);
275 
276  if(singleFitVertex.isValid()) {
277  if( pv.covariance(3,3) > 0. ) svhelper::updateVertexTime(singleFitVertex);
278  recoVertices.push_back(VTX(singleFitVertex));
279 
280 
281 #ifdef VTXDEBUG
282  const VTX & extVertex = recoVertices.back();
283  GlobalVector vtxDir = GlobalVector(extVertex.p4().X(),extVertex.p4().Y(),extVertex.p4().Z());
284 
285 
286  std::cout << " pointing : " << Geom::deltaR(extVertex.position() - pv.position(), vtxDir) << std::endl;
287 #endif
288  }
289  }
290  }
291  return recoVertices;
292 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
const double w
Definition: UKUtility.cc:23
void setBeamSpot(const reco::BeamSpot &beamSpot)
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:37
bool trackFilterArbitrator(const reco::TransientTrack &track) const
double timeExt() const
const MagneticField * field() const
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
const Point & position() const
position
Definition: Vertex.h:109
bool isFinite(T x)
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
double significance() const
Definition: Measurement1D.h:32
const Track & track() const
double dtErrorExt() const
double value() const
Definition: Measurement1D.h:28
void updateVertexTime(TransientVertex &vtx)
Definition: SVTimeHelpers.h:11
fixed size matrix
double cov33(const reco::Vertex &sv)
TrajectoryStateOnSurface impactPointState() const
bool isValid() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
float trackWeight(const reco::Vertex &sv, const reco::TransientTrack &track)
Global3DVector GlobalVector
Definition: GlobalVector.h:10

Member Data Documentation

template<class VTX >
edm::InputTag TrackVertexArbitration< VTX >::beamSpotCollection
private

Definition at line 77 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::distCut
private
template<class VTX >
double TrackVertexArbitration< VTX >::dLenFraction
private
template<class VTX >
double TrackVertexArbitration< VTX >::dRCut
private
template<class VTX >
double TrackVertexArbitration< VTX >::fitterRatio
private
template<class VTX >
double TrackVertexArbitration< VTX >::fitterSigmacut
private
template<class VTX >
double TrackVertexArbitration< VTX >::fitterTini
private
template<class VTX >
double TrackVertexArbitration< VTX >::maxTimeSignificance
private
template<class VTX >
edm::InputTag TrackVertexArbitration< VTX >::primaryVertexCollection
private

Definition at line 74 of file TrackVertexArbitration.h.

template<class VTX >
edm::InputTag TrackVertexArbitration< VTX >::secondaryVertexCollection
private

Definition at line 75 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::sigCut
private
template<class VTX >
double TrackVertexArbitration< VTX >::sign
private
template<class VTX >
edm::InputTag TrackVertexArbitration< VTX >::trackCollection
private

Definition at line 76 of file TrackVertexArbitration.h.

template<class VTX >
int TrackVertexArbitration< VTX >::trackMinLayers
private
template<class VTX >
int TrackVertexArbitration< VTX >::trackMinPixels
private
template<class VTX >
double TrackVertexArbitration< VTX >::trackMinPt
private