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);
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;
53 
54  std::auto_ptr<VertexReconstructor> vtxReco;
55  std::auto_ptr<TracksClusteringFromDisplacedSeed> clusterizer;
56 
57 };
58 
60  beamSpotCollection(params.getParameter<edm::InputTag>("beamSpot")),
61  primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
62  trackCollection(params.getParameter<edm::InputTag>("tracks")),
63  minHits(params.getParameter<unsigned int>("minHits")),
64  maxNTracks(params.getParameter<unsigned int>("maxNTracks")),
65  maxLIP(params.getParameter<double>("maximumLongitudinalImpactParameter")),
66  minPt(params.getParameter<double>("minPt")), //0.8
67  vertexMinAngleCosine(params.getParameter<double>("vertexMinAngleCosine")), //0.98
68  vertexMinDLen2DSig(params.getParameter<double>("vertexMinDLen2DSig")), //2.5
69  vertexMinDLenSig(params.getParameter<double>("vertexMinDLenSig")), //0.5
70  vtxReco(new ConfigurableVertexReconstructor(params.getParameter<edm::ParameterSet>("vertexReco"))),
71  clusterizer(new TracksClusteringFromDisplacedSeed(params.getParameter<edm::ParameterSet>("clusterizer")))
72 
73 {
74  produces<reco::VertexCollection>();
75  //produces<reco::VertexCollection>("multi");
76 }
77 
79 {
80  if (track->hitPattern().numberOfValidHits() < (int)minHits)
81 // if (track->hitPattern().trackerLayersWithMeasurement() < (int)minHits)
82  return false;
83  if (track->pt() < minPt )
84  return false;
85 
86  return true;
87 }
88 
90 {
91  using namespace reco;
92 
93  double sigmacut = 3.0;
94  double Tini = 256.;
95  double ratio = 0.25;
96  VertexDistance3D vdist;
97  VertexDistanceXY vdist2d;
98  MultiVertexFitter theMultiVertexFitter;
99  AdaptiveVertexFitter theAdaptiveFitter(
100  GeometricAnnealing(sigmacut, Tini, ratio),
105 
106 
108  event.getByLabel(beamSpotCollection, beamSpot);
109 
110  edm::Handle<VertexCollection> primaryVertices;
111  event.getByLabel(primaryVertexCollection, primaryVertices);
112 
114  event.getByLabel(trackCollection, tracks);
115 
117  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
118  trackBuilder);
119 
120 
121  std::auto_ptr<VertexCollection> recoVertices(new VertexCollection);
122  if(primaryVertices->size()!=0) {
123 
124  const reco::Vertex &pv = (*primaryVertices)[0];
125 
126  std::vector<TransientTrack> tts;
127  //Fill transient track vector
128  for(TrackCollection::const_iterator track = tracks->begin();
129  track != tracks->end(); ++track) {
130  TrackRef ref(tracks, track - tracks->begin());
131  if (!trackFilter(ref))
132  continue;
133  if( std::abs(ref->dz(pv.position())) > maxLIP)
134  continue;
135  TransientTrack tt = trackBuilder->build(ref);
136  tt.setBeamSpot(*beamSpot);
137  tts.push_back(tt);
138  }
139  std::vector<TracksClusteringFromDisplacedSeed::Cluster> clusters = clusterizer->clusters(pv,tts);
140 
141  //Create BS object from PV to feed in the AVR
143  for(unsigned int i = 0; i < 7; i++) {
144  for(unsigned int j = 0; j < 7; j++) {
145  if (i < 3 && j < 3)
146  cov(i, j) = pv.covariance(i, j);
147  else
148  cov(i, j) = 0.0;
149  }
150  }
151  BeamSpot bs(pv.position(), 0.0, 0.0, 0.0, 0.0, cov, BeamSpot::Unknown);
152 
153 
154  int i=0;
155 #ifdef VTXDEBUG
156 
157  std::cout << "CLUSTERS " << clusters.size() << std::endl;
158 #endif
159 
160  for(std::vector<TracksClusteringFromDisplacedSeed::Cluster>::iterator cluster = clusters.begin();
161  cluster != clusters.end(); ++cluster,++i)
162  {
163  if(cluster->tracks.size() == 0 || cluster->tracks.size() > maxNTracks )
164  continue;
165 
166  cluster->tracks.push_back(cluster->seedingTrack); //add the seed to the list of tracks to fit
167  std::vector<TransientVertex> vertices;
168  vertices = vtxReco->vertices(cluster->tracks, bs); // attempt with config given reconstructor
169  TransientVertex singleFitVertex;
170  singleFitVertex = theAdaptiveFitter.vertex(cluster->tracks,cluster->seedPoint); //attempt with direct fitting
171  if(singleFitVertex.isValid())
172  vertices.push_back(singleFitVertex);
173  for(std::vector<TransientVertex>::const_iterator v = vertices.begin();
174  v != vertices.end(); ++v) {
175 // if(v->degreesOfFreedom() > 0.2)
176  {
177  Measurement1D dlen= vdist.distance(pv,*v);
178  Measurement1D dlen2= vdist2d.distance(pv,*v);
179  reco::Vertex vv(*v);
180 #ifdef VTXDEBUG
181  std::cout << "V chi2/n: " << v->normalisedChiSquared() << " ndof: " <<v->degreesOfFreedom() ;
182  std::cout << " dlen: " << dlen.value() << " error: " << dlen.error() << " signif: " << dlen.significance();
183  std::cout << " dlen2: " << dlen2.value() << " error2: " << dlen2.error() << " signif2: " << dlen2.significance();
184  std::cout << " pos: " << vv.position() << " error: " <<vv.xError() << " " << vv.yError() << " " << vv.zError() << std::endl;
185 #endif
186  GlobalVector dir;
187  std::vector<reco::TransientTrack> ts = v->originalTracks();
188  for(std::vector<reco::TransientTrack>::const_iterator i = ts.begin();
189  i != ts.end(); ++i) {
190  reco::TrackRef t = i->trackBaseRef().castTo<reco::TrackRef>();
191  float w = v->trackWeight(*i);
192  if (w > 0.5) dir+=i->impactPointState().globalDirection();
193 #ifdef VTXDEBUG
194  std::cout << "\t[" << (*t).pt() << ": "
195  << (*t).eta() << ", "
196  << (*t).phi() << "], "
197  << w << std::endl;
198 #endif
199  }
200  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
201  GlobalPoint sv((*v).position().x(),(*v).position().y(),(*v).position().z());
202  float vscal = dir.unit().dot((sv-ppv).unit()) ;
203 // std::cout << "Vscal: " << vscal << std::endl;
204  if(dlen.significance() > vertexMinDLenSig && vscal > vertexMinAngleCosine && v->normalisedChiSquared() < 10 && dlen2.significance() > vertexMinDLen2DSig)
205  {
206  recoVertices->push_back(*v);
207 #ifdef VTXDEBUG
208 
209  std::cout << "ADDED" << std::endl;
210 #endif
211 
212  }
213  }
214  }
215  }
216 #ifdef VTXDEBUG
217 
218  std::cout << "Final put " << recoVertices->size() << std::endl;
219 #endif
220  }
221 
222  event.put(recoVertices);
223 
224 }
225 
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:32
int i
Definition: DBlmapReader.cc:9
void setBeamSpot(const reco::BeamSpot &beamSpot)
double zError() const
error on z
Definition: Vertex.h:105
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::InputTag primaryVertexCollection
double error() const
Definition: Measurement1D.h:30
#define abs(x)
Definition: mlp_lapack.h:159
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:110
virtual void produce(edm::Event &event, const edm::EventSetup &es)
const Point & position() const
position
Definition: Vertex.h:93
std::auto_ptr< TracksClusteringFromDisplacedSeed > clusterizer
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
int j
Definition: DBlmapReader.cc:9
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:101
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
T w() const
dbl *** dir
Definition: mlp_gen.cc:35
bool isValid() const
double yError() const
error on y
Definition: Vertex.h:103