CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
 
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 80 of file TrackVertexArbitratration.h.

References TrackVertexArbitration< VTX >::dRCut.

80  :
81  primaryVertexCollection (params.getParameter<edm::InputTag>("primaryVertices")),
82  secondaryVertexCollection (params.getParameter<edm::InputTag>("secondaryVertices")),
83  trackCollection (params.getParameter<edm::InputTag>("tracks")),
84  beamSpotCollection (params.getParameter<edm::InputTag>("beamSpot")),
85  dRCut (params.getParameter<double>("dRCut")),
86  distCut (params.getParameter<double>("distCut")),
87  sigCut (params.getParameter<double>("sigCut")),
88  dLenFraction (params.getParameter<double>("dLenFraction")),
89  fitterSigmacut (params.getParameter<double>("fitterSigmacut")),
90  fitterTini (params.getParameter<double>("fitterTini")),
91  fitterRatio (params.getParameter<double>("fitterRatio")),
92  trackMinLayers (params.getParameter<int32_t>("trackMinLayers")),
93  trackMinPt (params.getParameter<double>("trackMinPt")),
94  trackMinPixels (params.getParameter<int32_t>("trackMinPixels"))
95 {
96  dRCut*=dRCut;
97 }
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 99 of file TrackVertexArbitratration.h.

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

100 {
101  if(!track.isValid())
102  return false;
104  return false;
105  if (track.track().pt() < trackMinPt)
106  return false;
108  return false;
109 
110  return true;
111 }
bool isValid() const
Make the ReferenceCountingProxy method to check validity public.
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:516
double pt() const
track transverse momentum
Definition: TrackBase.h:608
const Track & track() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:437
int numberOfValidPixelHits() const
Definition: HitPattern.h:816
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 136 of file TrackVertexArbitratration.h.

References IPTools::absoluteImpactParameter3D(), RecoVertex::convertError(), RecoVertex::convertPos(), gather_cfg::cout, deltaR(), reco::deltaR2(), HLT_25ns14e33_v1_cff::distCut, HLT_25ns14e33_v1_cff::dLenFraction, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, HLT_25ns14e33_v1_cff::dRCut, reco::TransientTrack::field(), HLT_25ns14e33_v1_cff::fitterRatio, HLT_25ns14e33_v1_cff::fitterSigmacut, HLT_25ns14e33_v1_cff::fitterTini, reco::TransientTrack::impactPointState(), TransientVertex::isValid(), reco::Vertex::position(), dt_dqm_sourceclient_common_cff::reco, reco::Track::residualX(), reco::Track::residualY(), TrackCollections2monitor_cff::selectedTracks, reco::TransientTrack::setBeamSpot(), HLT_25ns14e33_v1_cff::sigCut, Measurement1D::significance(), reco::TransientTrack::track(), trackWeight(), groupFilesInBlocks::tt, Measurement1D::value(), AdaptiveVertexFitter::vertex(), and w.

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

Definition at line 67 of file TrackVertexArbitratration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::dLenFraction
private

Definition at line 69 of file TrackVertexArbitratration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::dRCut
private
template<class VTX >
double TrackVertexArbitration< VTX >::fitterRatio
private

Definition at line 72 of file TrackVertexArbitratration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::fitterSigmacut
private

Definition at line 70 of file TrackVertexArbitratration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::fitterTini
private

Definition at line 71 of file TrackVertexArbitratration.h.

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

Definition at line 68 of file TrackVertexArbitratration.h.

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

Definition at line 73 of file TrackVertexArbitratration.h.

template<class VTX >
int TrackVertexArbitration< VTX >::trackMinPixels
private

Definition at line 75 of file TrackVertexArbitratration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::trackMinPt
private

Definition at line 74 of file TrackVertexArbitratration.h.