CMS 3D CMS Logo

InclusiveVertexFinder.h
Go to the documentation of this file.
1 #ifndef InclusiveVertexFinder_h
2 #define InclusiveVertexFinder_h
3 #include <memory>
4 
10 
14 
21 
26 
32 
34 //#define VTXDEBUG 1
35 template <class InputContainer, class VTX>
37  public:
38  typedef std::vector<VTX> Product;
39  typedef typename InputContainer::value_type TRK;
41 
42  virtual void produce(edm::Event &event, const edm::EventSetup &es) override;
43 
44  private:
45  bool trackFilter(const reco::Track &track) const;
46  std::pair<std::vector<reco::TransientTrack>,GlobalPoint> nearTracks(const reco::TransientTrack &seed, const std::vector<reco::TransientTrack> & tracks, const reco::Vertex & primaryVertex) const;
47 
51  unsigned int minHits;
52  unsigned int maxNTracks;
53  double maxLIP;
54  double minPt;
59  double fitterTini;
60  double fitterRatio;
63  std::auto_ptr<VertexReconstructor> vtxReco;
64  std::auto_ptr<TracksClusteringFromDisplacedSeed> clusterizer;
65 
66 };
67 template <class InputContainer, class VTX>
69  minHits(params.getParameter<unsigned int>("minHits")),
70  maxNTracks(params.getParameter<unsigned int>("maxNTracks")),
71  maxLIP(params.getParameter<double>("maximumLongitudinalImpactParameter")),
72  minPt(params.getParameter<double>("minPt")), //0.8
73  vertexMinAngleCosine(params.getParameter<double>("vertexMinAngleCosine")), //0.98
74  vertexMinDLen2DSig(params.getParameter<double>("vertexMinDLen2DSig")), //2.5
75  vertexMinDLenSig(params.getParameter<double>("vertexMinDLenSig")), //0.5
76  fitterSigmacut(params.getParameter<double>("fitterSigmacut")),
77  fitterTini(params.getParameter<double>("fitterTini")),
78  fitterRatio(params.getParameter<double>("fitterRatio")),
79  useVertexFitter(params.getParameter<bool>("useDirectVertexFitter")),
80  useVertexReco(params.getParameter<bool>("useVertexReco")),
81  vtxReco(new ConfigurableVertexReconstructor(params.getParameter<edm::ParameterSet>("vertexReco"))),
82  clusterizer(new TracksClusteringFromDisplacedSeed(params.getParameter<edm::ParameterSet>("clusterizer")))
83 
84 {
85  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
86  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
87  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
88  produces<Product>();
89  //produces<reco::VertexCollection>("multi");
90 }
91 template <class InputContainer, class VTX>
93 {
94  if (track.hitPattern().numberOfValidHits() < (int)minHits)
95  return false;
96  if (track.pt() < minPt )
97  return false;
98 
99  return true;
100 }
101 
102 template <class InputContainer, class VTX>
104 {
105  using namespace reco;
106 
107  VertexDistance3D vdist;
108  VertexDistanceXY vdist2d;
109  MultiVertexFitter theMultiVertexFitter;
110  AdaptiveVertexFitter theAdaptiveFitter(
116 
117 
119  event.getByToken(token_beamSpot,beamSpot);
120 
122  event.getByToken(token_primaryVertex, primaryVertices);
123 
125  event.getByToken(token_tracks, tracks);
126 
128  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
129  trackBuilder);
130 
131 
132  auto recoVertices = std::make_unique<Product>();
133  if(primaryVertices->size()!=0) {
134 
135  const reco::Vertex &pv = (*primaryVertices)[0];
136  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
137 
138  std::vector<TransientTrack> tts;
139  //Fill transient track vector
140  for(typename InputContainer::const_iterator track = tracks->begin();
141  track != tracks->end(); ++track) {
142 //TransientTrack tt = trackBuilder->build(ref);
143 //TrackRef ref(tracks, track - tracks->begin());
144  TransientTrack tt(tthelpers::buildTT(tracks,trackBuilder, track - tracks->begin()));
145  if(!tt.isValid()) continue;
146  if (!trackFilter(tt.track()))
147  continue;
148  if( std::abs(tt.track().dz(pv.position())) > maxLIP)
149  continue;
150  tt.setBeamSpot(*beamSpot);
151  tts.push_back(tt);
152  }
153  std::vector<TracksClusteringFromDisplacedSeed::Cluster> clusters = clusterizer->clusters(pv,tts);
154 
155  //Create BS object from PV to feed in the AVR
157  for(unsigned int i = 0; i < 7; i++) {
158  for(unsigned int j = 0; j < 7; j++) {
159  if (i < 3 && j < 3)
160  cov(i, j) = pv.covariance(i, j);
161  else
162  cov(i, j) = 0.0;
163  }
164  }
165  BeamSpot bs(pv.position(), 0.0, 0.0, 0.0, 0.0, cov, BeamSpot::Unknown);
166 
167 
168  int i=0;
169 #ifdef VTXDEBUG
170 
171  std::cout << "CLUSTERS " << clusters.size() << std::endl;
172 #endif
173 
174  for(std::vector<TracksClusteringFromDisplacedSeed::Cluster>::iterator cluster = clusters.begin();
175  cluster != clusters.end(); ++cluster,++i)
176  {
177  if(cluster->tracks.size() < 2 || cluster->tracks.size() > maxNTracks )
178  continue;
179  std::vector<TransientVertex> vertices;
180  if(useVertexReco) {
181  vertices = vtxReco->vertices(cluster->tracks, bs); // attempt with config given reconstructor
182  }
183  TransientVertex singleFitVertex;
184  if(useVertexFitter) {
185  singleFitVertex = theAdaptiveFitter.vertex(cluster->tracks,cluster->seedPoint); //attempt with direct fitting
186  if(singleFitVertex.isValid())
187  vertices.push_back(singleFitVertex);
188  }
189  for(std::vector<TransientVertex>::const_iterator v = vertices.begin();
190  v != vertices.end(); ++v) {
191  Measurement1D dlen= vdist.distance(pv,*v);
192  Measurement1D dlen2= vdist2d.distance(pv,*v);
193 #ifdef VTXDEBUG
194  VTX vv(*v);
195  std::cout << "V chi2/n: " << v->normalisedChiSquared() << " ndof: " <<v->degreesOfFreedom() ;
196  std::cout << " dlen: " << dlen.value() << " error: " << dlen.error() << " signif: " << dlen.significance();
197  std::cout << " dlen2: " << dlen2.value() << " error2: " << dlen2.error() << " signif2: " << dlen2.significance();
198  std::cout << " pos: " << vv.position() << " error: " <<vv.xError() << " " << vv.yError() << " " << vv.zError() << std::endl;
199 #endif
200  GlobalVector dir;
201  std::vector<reco::TransientTrack> ts = v->originalTracks();
202  for(std::vector<reco::TransientTrack>::const_iterator i = ts.begin();
203  i != ts.end(); ++i) {
204  float w = v->trackWeight(*i);
205  if (w > 0.5) dir+=i->impactPointState().globalDirection();
206 #ifdef VTXDEBUG
207  std::cout << "\t[" << (*i).track().pt() << ": "
208  << (*i).track().eta() << ", "
209  << (*i).track().phi() << "], "
210  << w << std::endl;
211 #endif
212  }
213  GlobalPoint sv((*v).position().x(),(*v).position().y(),(*v).position().z());
214  float vscal = dir.unit().dot((sv-ppv).unit()) ;
215  // std::cout << "Vscal: " << vscal << std::endl;
216  if(dlen.significance() > vertexMinDLenSig && vscal > vertexMinAngleCosine && v->normalisedChiSquared() < 10 && dlen2.significance() > vertexMinDLen2DSig)
217  {
218  recoVertices->push_back(*v);
219 #ifdef VTXDEBUG
220  std::cout << "ADDED" << std::endl;
221 #endif
222  }
223 
224  }
225  }
226 #ifdef VTXDEBUG
227 
228  std::cout << "Final put " << recoVertices->size() << std::endl;
229 #endif
230  }
231 
232  event.put(std::move(recoVertices));
233 
234 }
235 #endif
virtual Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
T getParameter(std::string const &) const
std::auto_ptr< TracksClusteringFromDisplacedSeed > clusterizer
const double w
Definition: UKUtility.cc:23
void setBeamSpot(const reco::BeamSpot &beamSpot)
bool isValid() const
Make the ReferenceCountingProxy method to check validity public.
int numberOfValidHits() const
Definition: HitPattern.h:823
std::pair< std::vector< reco::TransientTrack >, GlobalPoint > nearTracks(const reco::TransientTrack &seed, const std::vector< reco::TransientTrack > &tracks, const reco::Vertex &primaryVertex) const
double error() const
Definition: Measurement1D.h:30
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
InputContainer::value_type TRK
const Point & position() const
position
Definition: Vertex.h:109
double pt() const
track transverse momentum
Definition: TrackBase.h:616
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< InputContainer > token_tracks
virtual Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const
bool trackFilter(const reco::Track &track) const
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:604
Vector3DBase unit() const
Definition: Vector3DBase.h:57
double significance() const
Definition: Measurement1D.h:32
virtual void produce(edm::Event &event, const edm::EventSetup &es) override
const Track & track() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
const T & get() const
Definition: EventSetup.h:56
double value() const
Definition: Measurement1D.h:28
fixed size matrix
HLT enums.
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
dbl *** dir
Definition: mlp_gen.cc:35
TemplatedInclusiveVertexFinder(const edm::ParameterSet &params)
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
std::auto_ptr< VertexReconstructor > vtxReco
reco::TransientTrack buildTT(edm::Handle< reco::TrackCollection > &tracks, edm::ESHandle< TransientTrackBuilder > &trackbuilder, unsigned int k)
Definition: TTHelpers.h:8
bool isValid() const
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1