CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
InclusiveVertexFinder.cc
Go to the documentation of this file.
1 #include <memory>
2 
8 
12 
19 
24 
30 
31 //#define VTXDEBUG 1
32 
34  public:
36 
37  virtual void produce(edm::Event &event, const edm::EventSetup &es) override;
38 
39  private:
40  bool trackFilter(const reco::TrackRef &track) const;
41  std::pair<std::vector<reco::TransientTrack>,GlobalPoint> nearTracks(const reco::TransientTrack &seed, const std::vector<reco::TransientTrack> & tracks, const reco::Vertex & primaryVertex) const;
42 
46  unsigned int minHits;
47  unsigned int maxNTracks;
48  double maxLIP;
49  double minPt;
54  double fitterTini;
55  double fitterRatio;
58  std::auto_ptr<VertexReconstructor> vtxReco;
59  std::auto_ptr<TracksClusteringFromDisplacedSeed> clusterizer;
60 
61 };
62 
64  minHits(params.getParameter<unsigned int>("minHits")),
65  maxNTracks(params.getParameter<unsigned int>("maxNTracks")),
66  maxLIP(params.getParameter<double>("maximumLongitudinalImpactParameter")),
67  minPt(params.getParameter<double>("minPt")), //0.8
68  vertexMinAngleCosine(params.getParameter<double>("vertexMinAngleCosine")), //0.98
69  vertexMinDLen2DSig(params.getParameter<double>("vertexMinDLen2DSig")), //2.5
70  vertexMinDLenSig(params.getParameter<double>("vertexMinDLenSig")), //0.5
71  fitterSigmacut(params.getParameter<double>("fitterSigmacut")),
72  fitterTini(params.getParameter<double>("fitterTini")),
73  fitterRatio(params.getParameter<double>("fitterRatio")),
74  useVertexFitter(params.getParameter<bool>("useDirectVertexFitter")),
75  useVertexReco(params.getParameter<bool>("useVertexReco")),
76  vtxReco(new ConfigurableVertexReconstructor(params.getParameter<edm::ParameterSet>("vertexReco"))),
77  clusterizer(new TracksClusteringFromDisplacedSeed(params.getParameter<edm::ParameterSet>("clusterizer")))
78 
79 {
80  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
81  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
82  token_tracks = consumes<reco::TrackCollection>(params.getParameter<edm::InputTag>("tracks"));
83  produces<reco::VertexCollection>();
84  //produces<reco::VertexCollection>("multi");
85 }
86 
88 {
89  if (track->hitPattern().numberOfValidHits() < (int)minHits)
90 // if (track->hitPattern().trackerLayersWithMeasurement() < (int)minHits)
91  return false;
92  if (track->pt() < minPt )
93  return false;
94 
95  return true;
96 }
97 
99 {
100  using namespace reco;
101 
102  VertexDistance3D vdist;
103  VertexDistanceXY vdist2d;
104  MultiVertexFitter theMultiVertexFitter;
105  AdaptiveVertexFitter theAdaptiveFitter(
111 
112 
114  event.getByToken(token_beamSpot,beamSpot);
115 
116  edm::Handle<VertexCollection> primaryVertices;
117  event.getByToken(token_primaryVertex, primaryVertices);
118 
120  event.getByToken(token_tracks, tracks);
121 
123  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
124  trackBuilder);
125 
126 
127  std::auto_ptr<VertexCollection> recoVertices(new VertexCollection);
128  if(primaryVertices->size()!=0) {
129 
130  const reco::Vertex &pv = (*primaryVertices)[0];
131 
132  std::vector<TransientTrack> tts;
133  //Fill transient track vector
134  for(TrackCollection::const_iterator track = tracks->begin();
135  track != tracks->end(); ++track) {
136  TrackRef ref(tracks, track - tracks->begin());
137  if (!trackFilter(ref))
138  continue;
139  if( std::abs(ref->dz(pv.position())) > maxLIP)
140  continue;
141  TransientTrack tt = trackBuilder->build(ref);
142  tt.setBeamSpot(*beamSpot);
143  tts.push_back(tt);
144  }
145  std::vector<TracksClusteringFromDisplacedSeed::Cluster> clusters = clusterizer->clusters(pv,tts);
146 
147  //Create BS object from PV to feed in the AVR
149  for(unsigned int i = 0; i < 7; i++) {
150  for(unsigned int j = 0; j < 7; j++) {
151  if (i < 3 && j < 3)
152  cov(i, j) = pv.covariance(i, j);
153  else
154  cov(i, j) = 0.0;
155  }
156  }
157  BeamSpot bs(pv.position(), 0.0, 0.0, 0.0, 0.0, cov, BeamSpot::Unknown);
158 
159 
160  int i=0;
161 #ifdef VTXDEBUG
162 
163  std::cout << "CLUSTERS " << clusters.size() << std::endl;
164 #endif
165 
166  for(std::vector<TracksClusteringFromDisplacedSeed::Cluster>::iterator cluster = clusters.begin();
167  cluster != clusters.end(); ++cluster,++i)
168  {
169  if(cluster->tracks.size() < 2 || cluster->tracks.size() > maxNTracks )
170  continue;
171  std::vector<TransientVertex> vertices;
172  if(useVertexReco) {
173  vertices = vtxReco->vertices(cluster->tracks, bs); // attempt with config given reconstructor
174  }
175  TransientVertex singleFitVertex;
176  if(useVertexFitter) {
177  singleFitVertex = theAdaptiveFitter.vertex(cluster->tracks,cluster->seedPoint); //attempt with direct fitting
178  if(singleFitVertex.isValid())
179  vertices.push_back(singleFitVertex);
180  }
181  for(std::vector<TransientVertex>::const_iterator v = vertices.begin();
182  v != vertices.end(); ++v) {
183 // if(v->degreesOfFreedom() > 0.2)
184  {
185  Measurement1D dlen= vdist.distance(pv,*v);
186  Measurement1D dlen2= vdist2d.distance(pv,*v);
187  reco::Vertex vv(*v);
188 #ifdef VTXDEBUG
189  std::cout << "V chi2/n: " << v->normalisedChiSquared() << " ndof: " <<v->degreesOfFreedom() ;
190  std::cout << " dlen: " << dlen.value() << " error: " << dlen.error() << " signif: " << dlen.significance();
191  std::cout << " dlen2: " << dlen2.value() << " error2: " << dlen2.error() << " signif2: " << dlen2.significance();
192  std::cout << " pos: " << vv.position() << " error: " <<vv.xError() << " " << vv.yError() << " " << vv.zError() << std::endl;
193 #endif
194  GlobalVector dir;
195  std::vector<reco::TransientTrack> ts = v->originalTracks();
196  for(std::vector<reco::TransientTrack>::const_iterator i = ts.begin();
197  i != ts.end(); ++i) {
198  reco::TrackRef t = i->trackBaseRef().castTo<reco::TrackRef>();
199  float w = v->trackWeight(*i);
200  if (w > 0.5) dir+=i->impactPointState().globalDirection();
201 #ifdef VTXDEBUG
202  std::cout << "\t[" << (*t).pt() << ": "
203  << (*t).eta() << ", "
204  << (*t).phi() << "], "
205  << w << std::endl;
206 #endif
207  }
208  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
209  GlobalPoint sv((*v).position().x(),(*v).position().y(),(*v).position().z());
210  float vscal = dir.unit().dot((sv-ppv).unit()) ;
211 // std::cout << "Vscal: " << vscal << std::endl;
212  if(dlen.significance() > vertexMinDLenSig && vscal > vertexMinAngleCosine && v->normalisedChiSquared() < 10 && dlen2.significance() > vertexMinDLen2DSig)
213  {
214  recoVertices->push_back(*v);
215 #ifdef VTXDEBUG
216 
217  std::cout << "ADDED" << std::endl;
218 #endif
219 
220  }
221  }
222  }
223  }
224 #ifdef VTXDEBUG
225 
226  std::cout << "Final put " << recoVertices->size() << std::endl;
227 #endif
228  }
229 
230  event.put(recoVertices);
231 
232 }
233 
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
const double w
Definition: UKUtility.cc:23
void setBeamSpot(const reco::BeamSpot &beamSpot)
double zError() const
error on z
Definition: Vertex.h:118
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double error() const
Definition: Measurement1D.h:30
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
std::pair< std::vector< reco::TransientTrack >, GlobalPoint > nearTracks(const reco::TransientTrack &seed, const std::vector< reco::TransientTrack > &tracks, const reco::Vertex &primaryVertex) const
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:123
const Point & position() const
position
Definition: Vertex.h:106
std::auto_ptr< TracksClusteringFromDisplacedSeed > clusterizer
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
virtual void produce(edm::Event &event, const edm::EventSetup &es) override
bool trackFilter(const reco::TrackRef &track) const
string unit
Definition: csvLumiCalc.py:46
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
std::auto_ptr< VertexReconstructor > vtxReco
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::TrackCollection > token_tracks
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
InclusiveVertexFinder(const edm::ParameterSet &params)
Vector3DBase unit() const
Definition: Vector3DBase.h:57
double significance() const
Definition: Measurement1D.h:32
double xError() const
error on x
Definition: Vertex.h:114
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:55
double value() const
Definition: Measurement1D.h:28
tuple cout
Definition: gather_cfg.py:121
dbl *** dir
Definition: mlp_gen.cc:35
bool isValid() const
double yError() const
error on y
Definition: Vertex.h:116
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex