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 
36 
37 #include <type_traits>
38 
39 //#define VTXDEBUG 1
40 template <class InputContainer, class VTX>
42  public:
43  typedef std::vector<VTX> Product;
44  typedef typename InputContainer::value_type TRK;
46 
49  pdesc.add<edm::InputTag>("beamSpot",edm::InputTag("offlineBeamSpot"));
50  pdesc.add<edm::InputTag>("primaryVertices",edm::InputTag("offlinePrimaryVertices"));
52  pdesc.add<edm::InputTag>("tracks",edm::InputTag("generalTracks"));
53  pdesc.add<unsigned int>("minHits",8);
55  pdesc.add<edm::InputTag>("tracks",edm::InputTag("particleFlow"));
56  pdesc.add<unsigned int>("minHits",0);
57  } else {
58  pdesc.add<edm::InputTag>("tracks",edm::InputTag("generalTracks"));
59  }
60 
61  pdesc.add<double>("maximumLongitudinalImpactParameter",0.3);
62  pdesc.add<double>("maximumTimeSignificance",3.0);
63  pdesc.add<double>("minPt",0.8);
64  pdesc.add<unsigned int>("maxNTracks",30);
65  //clusterizer pset
67  clusterizer.add<double>("seedMax3DIPSignificance",9999.0);
68  clusterizer.add<double>("seedMax3DIPValue",9999.0);
69  clusterizer.add<double>("seedMin3DIPSignificance",1.2);
70  clusterizer.add<double>("seedMin3DIPValue",0.005);
71  clusterizer.add<double>("clusterMaxDistance",0.05);
72  clusterizer.add<double>("clusterMaxSignificance",4.5);
73  clusterizer.add<double>("distanceRatio",20.0);
74  clusterizer.add<double>("clusterMinAngleCosine",0.5);
75  clusterizer.add<double>("maxTimeSignificance",3.5);
76  pdesc.add<edm::ParameterSetDescription>("clusterizer", clusterizer);
77  // vertex and fitter config
78  pdesc.add<double>("vertexMinAngleCosine",0.95);
79  pdesc.add<double>("vertexMinDLen2DSig",2.5);
80  pdesc.add<double>("vertexMinDLenSig",0.5);
81  pdesc.add<double>("fitterSigmacut",3.0);
82  pdesc.add<double>("fitterTini",256.0);
83  pdesc.add<double>("fitterRatio",0.25);
84  pdesc.add<bool>("useDirectVertexFitter",true);
85  pdesc.add<bool>("useVertexReco",true);
86  // vertexReco pset
88  vertexReco.add<std::string>("finder", std::string("avr"));
89  vertexReco.add<double>("primcut",1.0);
90  vertexReco.add<double>("seccut",3.0);
91  vertexReco.add<bool>("smoothing",true);
92  pdesc.add<edm::ParameterSetDescription>("vertexReco", vertexReco);
94  cdesc.add("inclusiveVertexFinderDefault",pdesc);
96  cdesc.add("inclusiveCandidateVertexFinderDefault",pdesc);
97  } else {
98  cdesc.addDefault(pdesc);
99  }
100  }
101 
102  virtual void produce(edm::Event &event, const edm::EventSetup &es) override;
103 
104  private:
105  bool trackFilter(const reco::Track &track) const;
106  std::pair<std::vector<reco::TransientTrack>,GlobalPoint> nearTracks(const reco::TransientTrack &seed, const std::vector<reco::TransientTrack> & tracks, const reco::Vertex & primaryVertex) const;
107 
111  unsigned int minHits;
112  unsigned int maxNTracks;
113  double maxLIP;
114  double maxTimeSig;
115  double minPt;
120  double fitterTini;
121  double fitterRatio;
124  std::auto_ptr<VertexReconstructor> vtxReco;
125  std::auto_ptr<TracksClusteringFromDisplacedSeed> clusterizer;
126 
127 };
128 template <class InputContainer, class VTX>
130  minHits(params.getParameter<unsigned int>("minHits")),
131  maxNTracks(params.getParameter<unsigned int>("maxNTracks")),
132  maxLIP(params.getParameter<double>("maximumLongitudinalImpactParameter")),
133  maxTimeSig(params.getParameter<double>("maximumTimeSignificance")),
134  minPt(params.getParameter<double>("minPt")), //0.8
135  vertexMinAngleCosine(params.getParameter<double>("vertexMinAngleCosine")), //0.98
136  vertexMinDLen2DSig(params.getParameter<double>("vertexMinDLen2DSig")), //2.5
137  vertexMinDLenSig(params.getParameter<double>("vertexMinDLenSig")), //0.5
138  fitterSigmacut(params.getParameter<double>("fitterSigmacut")),
139  fitterTini(params.getParameter<double>("fitterTini")),
140  fitterRatio(params.getParameter<double>("fitterRatio")),
141  useVertexFitter(params.getParameter<bool>("useDirectVertexFitter")),
142  useVertexReco(params.getParameter<bool>("useVertexReco")),
143  vtxReco(new ConfigurableVertexReconstructor(params.getParameter<edm::ParameterSet>("vertexReco"))),
144  clusterizer(new TracksClusteringFromDisplacedSeed(params.getParameter<edm::ParameterSet>("clusterizer")))
145 
146 {
147  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
148  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
149  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
150  produces<Product>();
151  //produces<reco::VertexCollection>("multi");
152 }
153 template <class InputContainer, class VTX>
155 {
156  if (track.hitPattern().numberOfValidHits() < (int)minHits)
157  return false;
158  if (track.pt() < minPt )
159  return false;
160 
161  return true;
162 }
163 
164 template <class InputContainer, class VTX>
166 {
167  using namespace reco;
168 
169  VertexDistance3D vdist;
170  VertexDistanceXY vdist2d;
171  MultiVertexFitter theMultiVertexFitter;
172  AdaptiveVertexFitter theAdaptiveFitter(
178 
179 
181  event.getByToken(token_beamSpot,beamSpot);
182 
184  event.getByToken(token_primaryVertex, primaryVertices);
185 
187  event.getByToken(token_tracks, tracks);
188 
190  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
191  trackBuilder);
192 
193 
194  auto recoVertices = std::make_unique<Product>();
195  if(primaryVertices->size()!=0) {
196 
197  const reco::Vertex &pv = (*primaryVertices)[0];
198  GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
199 
200  std::vector<TransientTrack> tts;
201  //Fill transient track vector
202  for(typename InputContainer::const_iterator track = tracks->begin();
203  track != tracks->end(); ++track) {
204 //TransientTrack tt = trackBuilder->build(ref);
205 //TrackRef ref(tracks, track - tracks->begin());
206  TransientTrack tt(tthelpers::buildTT(tracks,trackBuilder, track - tracks->begin()));
207  if(!tt.isValid()) continue;
208  if (!trackFilter(tt.track()))
209  continue;
210  if( std::abs(tt.track().dz(pv.position())) > maxLIP)
211  continue;
212  if( edm::isFinite(tt.timeExt()) && pv.covariance(3,3) > 0. ) { // only apply if time available
213  auto tError = std::sqrt( std::pow(tt.dtErrorExt(),2) + pv.covariance(3,3) );
214  auto dtSig = std::abs(tt.timeExt() - pv.t())/tError;
215  if( dtSig > maxTimeSig ) continue;
216  }
217  tt.setBeamSpot(*beamSpot);
218  tts.push_back(tt);
219  }
220  std::vector<TracksClusteringFromDisplacedSeed::Cluster> clusters = clusterizer->clusters(pv,tts);
221 
222  //Create BS object from PV to feed in the AVR
224  for(unsigned int i = 0; i < 7; i++) {
225  for(unsigned int j = 0; j < 7; j++) {
226  if (i < 3 && j < 3)
227  cov(i, j) = pv.covariance(i, j);
228  else
229  cov(i, j) = 0.0;
230  }
231  }
232  BeamSpot bs(pv.position(), 0.0, 0.0, 0.0, 0.0, cov, BeamSpot::Unknown);
233 
234 
235  int i=0;
236 #ifdef VTXDEBUG
237 
238  std::cout << "CLUSTERS " << clusters.size() << std::endl;
239 #endif
240 
241  for(std::vector<TracksClusteringFromDisplacedSeed::Cluster>::iterator cluster = clusters.begin();
242  cluster != clusters.end(); ++cluster,++i)
243  {
244  if(cluster->tracks.size() < 2 || cluster->tracks.size() > maxNTracks )
245  continue;
246  std::vector<TransientVertex> vertices;
247  if(useVertexReco) {
248  vertices = vtxReco->vertices(cluster->tracks, bs); // attempt with config given reconstructor
249  }
250  TransientVertex singleFitVertex;
251  if(useVertexFitter) {
252  singleFitVertex = theAdaptiveFitter.vertex(cluster->tracks,cluster->seedPoint); //attempt with direct fitting
253  if(singleFitVertex.isValid())
254  vertices.push_back(singleFitVertex);
255  }
256 
257  // for each transient vertex state determine if a time can be measured and fill covariance
258  if( pv.covariance(3,3) > 0. ) {
259  for(auto& vtx : vertices) {
261  }
262  }
263 
264  for(std::vector<TransientVertex>::const_iterator v = vertices.begin();
265  v != vertices.end(); ++v) {
266  Measurement1D dlen= vdist.distance(pv,*v);
267  Measurement1D dlen2= vdist2d.distance(pv,*v);
268 #ifdef VTXDEBUG
269  VTX vv(*v);
270  std::cout << "V chi2/n: " << v->normalisedChiSquared() << " ndof: " <<v->degreesOfFreedom() ;
271  std::cout << " dlen: " << dlen.value() << " error: " << dlen.error() << " signif: " << dlen.significance();
272  std::cout << " dlen2: " << dlen2.value() << " error2: " << dlen2.error() << " signif2: " << dlen2.significance();
273  std::cout << " pos: " << vv.position() << " error: " <<vv.xError() << " " << vv.yError() << " " << vv.zError() << std::endl;
274  std::cout << " time: " << vv.time() << " error: " << vv.tError() << std::endl;
275 #endif
276  GlobalVector dir;
277  std::vector<reco::TransientTrack> ts = v->originalTracks();
278  for(std::vector<reco::TransientTrack>::const_iterator i = ts.begin();
279  i != ts.end(); ++i) {
280  float w = v->trackWeight(*i);
281  if (w > 0.5) dir+=i->impactPointState().globalDirection();
282 #ifdef VTXDEBUG
283  std::cout << "\t[" << (*i).track().pt() << ": "
284  << (*i).track().eta() << ", "
285  << (*i).track().phi() << "], "
286  << w << std::endl;
287 #endif
288  }
289  GlobalPoint sv((*v).position().x(),(*v).position().y(),(*v).position().z());
290  float vscal = dir.unit().dot((sv-ppv).unit());
291  if(dlen.significance() > vertexMinDLenSig &&
292  ( (vertexMinAngleCosine > 0) ? (vscal > vertexMinAngleCosine) : (vscal < vertexMinAngleCosine) )
293  && v->normalisedChiSquared() < 10 && dlen2.significance() > vertexMinDLen2DSig)
294  {
295  recoVertices->push_back(*v);
296 
297 #ifdef VTXDEBUG
298  std::cout << "ADDED" << std::endl;
299 #endif
300  }
301 
302  }
303  }
304 #ifdef VTXDEBUG
305 
306  std::cout << "Final put " << recoVertices->size() << std::endl;
307 #endif
308  }
309 
310  event.put(std::move(recoVertices));
311 
312 }
313 #endif
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:824
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
double timeExt() const
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
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
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
bool isFinite(T x)
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
void addDefault(ParameterSetDescription const &psetDescription)
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:621
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< InputContainer > token_tracks
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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:609
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
primaryVertices
Definition: jets_cff.py:130
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:446
double dtErrorExt() const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double value() const
Definition: Measurement1D.h:28
void updateVertexTime(TransientVertex &vtx)
Definition: SVTimeHelpers.h:11
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:63
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:10
static void fillDescriptions(edm::ConfigurationDescriptions &cdesc)
bool isValid() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1
double t() const
t coordinate
Definition: Vertex.h:117