CMS 3D CMS Logo

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

#include <TrackVertexArbitratration.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
 
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 46 of file TrackVertexArbitratration.h.

Constructor & Destructor Documentation

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

Definition at line 81 of file TrackVertexArbitratration.h.

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

81  :
82  primaryVertexCollection (params.getParameter<edm::InputTag>("primaryVertices")),
83  secondaryVertexCollection (params.getParameter<edm::InputTag>("secondaryVertices")),
84  trackCollection (params.getParameter<edm::InputTag>("tracks")),
85  beamSpotCollection (params.getParameter<edm::InputTag>("beamSpot")),
86  dRCut (params.getParameter<double>("dRCut")),
87  distCut (params.getParameter<double>("distCut")),
88  sigCut (params.getParameter<double>("sigCut")),
89  dLenFraction (params.getParameter<double>("dLenFraction")),
90  fitterSigmacut (params.getParameter<double>("fitterSigmacut")),
91  fitterTini (params.getParameter<double>("fitterTini")),
92  fitterRatio (params.getParameter<double>("fitterRatio")),
93  trackMinLayers (params.getParameter<int32_t>("trackMinLayers")),
94  trackMinPt (params.getParameter<double>("trackMinPt")),
95  trackMinPixels (params.getParameter<int32_t>("trackMinPixels"))
96 {
97  sign = 1.0;
98  if(dRCut < 0){sign = -1.0;}
99  dRCut*=dRCut;
100 }
T getParameter(std::string const &) const

Member Function Documentation

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

Definition at line 102 of file TrackVertexArbitratration.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().

103 {
104  if(!track.isValid())
105  return false;
107  return false;
108  if (track.track().pt() < trackMinPt)
109  return false;
111  return false;
112 
113  return true;
114 }
bool isValid() const
Make the ReferenceCountingProxy method to check validity public.
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:518
double pt() const
track transverse momentum
Definition: TrackBase.h:621
const Track & track() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:446
int numberOfValidPixelHits() const
Definition: HitPattern.h:838
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 139 of file TrackVertexArbitratration.h.

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

144 {
145  using namespace reco;
146 
147  //std::cout << "PV: " << pv.position() << std::endl;
148  VertexDistance3D dist;
149  AdaptiveVertexFitter theAdaptiveFitter(
155 
156 
157 
158  std::vector<VTX> recoVertices;
159  VertexDistance3D vdist;
160  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
161 
162  std::map<unsigned int, Measurement1D> cachedIP;
163  for(typename std::vector<VTX>::const_iterator sv = secondaryVertices.begin();
164  sv != secondaryVertices.end(); ++sv) {
165 
166  GlobalPoint ssv(sv->position().x(),sv->position().y(),sv->position().z());
167  GlobalVector flightDir = ssv-ppv;
168 // std::cout << "Vertex : " << sv-secondaryVertices->begin() << " " << sv->position() << std::endl;
170  std::vector<reco::TransientTrack> selTracks;
171  for(unsigned int itrack = 0; itrack < selectedTracks.size(); itrack++){
172  TransientTrack & tt = (selectedTracks)[itrack];
173  if (!trackFilterArbitrator(tt)) continue;
174  tt.setBeamSpot(*beamSpot);
175  float w = trackWeight(*sv,tt);
176  Measurement1D ipv;
177  if(cachedIP.find(itrack)!=cachedIP.end()) {
178  ipv=cachedIP[itrack];
179  } else {
180  std::pair<bool,Measurement1D> ipvp = IPTools::absoluteImpactParameter3D(tt,pv);
181  cachedIP[itrack]=ipvp.second;
182  ipv=ipvp.second;
183  }
184 
185  AnalyticalImpactPointExtrapolator extrapolator(tt.field());
186  TrajectoryStateOnSurface tsos = extrapolator.extrapolate(tt.impactPointState(), RecoVertex::convertPos(sv->position()));
187  if(! tsos.isValid()) continue;
188  GlobalPoint refPoint = tsos.globalPosition();
189  GlobalError refPointErr = tsos.cartesianError().position();
190  Measurement1D isv = dist.distance(VertexState(RecoVertex::convertPos(sv->position()),RecoVertex::convertError(sv->error())),VertexState(refPoint, refPointErr));
191 
192  float dR = Geom::deltaR2(( (sign > 0) ? flightDir : -flightDir),tt.track()); //.eta(), flightDir.phi(), tt.track().eta(), tt.track().phi());
193 
194  if( w > 0 || ( isv.significance() < sigCut && isv.value() < distCut && isv.value() < dlen.value()*dLenFraction ) )
195  {
196 
197  if(( isv.value() < ipv.value() ) && isv.value() < distCut && isv.value() < dlen.value()*dLenFraction
198  && dR < dRCut )
199  {
200 #ifdef VTXDEBUG
201  if(w > 0.5) std::cout << " = ";
202  else std::cout << " + ";
203 #endif
204  selTracks.push_back(tt);
205  } else
206  {
207 #ifdef VTXDEBUG
208  if(w > 0.5 && isv.value() > ipv.value() ) std::cout << " - ";
209  else std::cout << " ";
210 #endif
211  //add also the tracks used in previous fitting that are still closer to Sv than Pv
212  if(w > 0.5 && isv.value() <= ipv.value() && dR < dRCut) {
213  selTracks.push_back(tt);
214 #ifdef VTXDEBUG
215  std::cout << " = ";
216 #endif
217  }
218  if(w > 0.5 && isv.value() <= ipv.value() && dR >= dRCut) {
219 #ifdef VTXDEBUG
220  std::cout << " - ";
221 #endif
222 
223  }
224 
225 
226  }
227 #ifdef VTXDEBUG
228 
229  std::cout << "t : " << itrack << " ref " << ref.key() << " w: " << w
230  << " svip: " << isv.significance() << " " << isv.value()
231  << " pvip: " << ipv.significance() << " " << ipv.value() << " res " << tt.track().residualX(0) << "," << tt.track().residualY(0)
232 // << " tpvip: " << itpv.second.significance() << " " << itpv.second.value() << " dr: " << dR << std::endl;
233 #endif
234 
235  }
236  else
237  {
238 #ifdef VTXDEBUG
239 
240  std::cout << " . t : " << itrack << " ref " << ref.key()<< " w: " << w
241  << " svip: " << isv.second.significance() << " " << isv.second.value()
242  << " pvip: " << ipv.significance() << " " << ipv.value() << " dr: " << dR << std::endl;
243 #endif
244 
245  }
246  }
247 
248  if(selTracks.size() >= 2)
249  {
250  TransientVertex singleFitVertex;
251  singleFitVertex = theAdaptiveFitter.vertex(selTracks,ssv);
252  if(singleFitVertex.isValid()) { recoVertices.push_back(VTX(singleFitVertex));
253 
254 #ifdef VTXDEBUG
255  const VTX & extVertex = recoVertices.back();
256  GlobalVector vtxDir = GlobalVector(extVertex.p4().X(),extVertex.p4().Y(),extVertex.p4().Z());
257 
258 
259  std::cout << " pointing : " << Geom::deltaR(extVertex.position() - pv.position(), vtxDir) << std::endl;
260 #endif
261  }
262  }
263  }
264  return recoVertices;
265 }
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
const MagneticField * field() const
const Point & position() const
position
Definition: Vertex.h:109
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double significance() const
Definition: Measurement1D.h:32
const Track & track() const
float trackWeight(const reco::Vertex &sv, const reco::TransientTrack &track)
double value() const
Definition: Measurement1D.h:28
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
fixed size matrix
TrajectoryStateOnSurface impactPointState() const
bool isValid() const
Global3DVector GlobalVector
Definition: GlobalVector.h:10

Member Data Documentation

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

Definition at line 65 of file TrackVertexArbitratration.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 >
edm::InputTag TrackVertexArbitration< VTX >::primaryVertexCollection
private

Definition at line 62 of file TrackVertexArbitratration.h.

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

Definition at line 63 of file TrackVertexArbitratration.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 64 of file TrackVertexArbitratration.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