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  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,
107  const std::vector<reco::TransientTrack> &tracks,
108  const reco::Vertex &primaryVertex) const;
109 
113  unsigned int minHits;
114  unsigned int maxNTracks;
115  double maxLIP;
116  double maxTimeSig;
117  double minPt;
122  double fitterTini;
123  double fitterRatio;
126  std::unique_ptr<VertexReconstructor> vtxReco;
127  std::unique_ptr<TracksClusteringFromDisplacedSeed> clusterizer;
128 };
129 template <class InputContainer, class VTX>
131  : minHits(params.getParameter<unsigned int>("minHits")),
132  maxNTracks(params.getParameter<unsigned int>("maxNTracks")),
133  maxLIP(params.getParameter<double>("maximumLongitudinalImpactParameter")),
134  maxTimeSig(params.getParameter<double>("maximumTimeSignificance")),
135  minPt(params.getParameter<double>("minPt")), //0.8
136  vertexMinAngleCosine(params.getParameter<double>("vertexMinAngleCosine")), //0.98
137  vertexMinDLen2DSig(params.getParameter<double>("vertexMinDLen2DSig")), //2.5
138  vertexMinDLenSig(params.getParameter<double>("vertexMinDLenSig")), //0.5
139  fitterSigmacut(params.getParameter<double>("fitterSigmacut")),
140  fitterTini(params.getParameter<double>("fitterTini")),
141  fitterRatio(params.getParameter<double>("fitterRatio")),
142  useVertexFitter(params.getParameter<bool>("useDirectVertexFitter")),
143  useVertexReco(params.getParameter<bool>("useVertexReco")),
144  vtxReco(new ConfigurableVertexReconstructor(params.getParameter<edm::ParameterSet>("vertexReco"))),
145  clusterizer(new TracksClusteringFromDisplacedSeed(params.getParameter<edm::ParameterSet>("clusterizer")))
146 
147 {
148  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
149  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
150  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
151  produces<Product>();
152  //produces<reco::VertexCollection>("multi");
153 }
154 template <class InputContainer, class VTX>
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  using namespace reco;
167 
168  VertexDistance3D vdist;
169  VertexDistanceXY vdist2d;
170  MultiVertexFitter theMultiVertexFitter;
176 
178  event.getByToken(token_beamSpot, beamSpot);
179 
181  event.getByToken(token_primaryVertex, primaryVertices);
182 
184  event.getByToken(token_tracks, tracks);
185 
187  es.get<TransientTrackRecord>().get("TransientTrackBuilder", trackBuilder);
188 
189  auto recoVertices = std::make_unique<Product>();
190  if (!primaryVertices->empty()) {
191  const reco::Vertex &pv = (*primaryVertices)[0];
192  GlobalPoint ppv(pv.position().x(), pv.position().y(), pv.position().z());
193 
194  std::vector<TransientTrack> tts;
195  //Fill transient track vector
196  for (typename InputContainer::const_iterator track = tracks->begin(); track != tracks->end(); ++track) {
197  //TransientTrack tt = trackBuilder->build(ref);
198  //TrackRef ref(tracks, track - tracks->begin());
199  TransientTrack tt(tthelpers::buildTT(tracks, trackBuilder, track - tracks->begin()));
200  if (!tt.isValid())
201  continue;
202  if (!trackFilter(tt.track()))
203  continue;
204  if (std::abs(tt.track().dz(pv.position())) > maxLIP)
205  continue;
206  if (edm::isFinite(tt.timeExt()) && pv.covariance(3, 3) > 0.) { // only apply if time available
207  auto tError = std::sqrt(std::pow(tt.dtErrorExt(), 2) + pv.covariance(3, 3));
208  auto dtSig = std::abs(tt.timeExt() - pv.t()) / tError;
209  if (dtSig > maxTimeSig)
210  continue;
211  }
212  tt.setBeamSpot(*beamSpot);
213  tts.push_back(tt);
214  }
215  std::vector<TracksClusteringFromDisplacedSeed::Cluster> clusters = clusterizer->clusters(pv, tts);
216 
217  //Create BS object from PV to feed in the AVR
219  for (unsigned int i = 0; i < 7; i++) {
220  for (unsigned int j = 0; j < 7; j++) {
221  if (i < 3 && j < 3)
222  cov(i, j) = pv.covariance(i, j);
223  else
224  cov(i, j) = 0.0;
225  }
226  }
227  BeamSpot bs(pv.position(), 0.0, 0.0, 0.0, 0.0, cov, BeamSpot::Unknown);
228 
229  int i = 0;
230 #ifdef VTXDEBUG
231 
232  std::cout << "CLUSTERS " << clusters.size() << std::endl;
233 #endif
234 
235  for (std::vector<TracksClusteringFromDisplacedSeed::Cluster>::iterator cluster = clusters.begin();
236  cluster != clusters.end();
237  ++cluster, ++i) {
238  if (cluster->tracks.size() < 2 || cluster->tracks.size() > maxNTracks)
239  continue;
240  std::vector<TransientVertex> vertices;
241  if (useVertexReco) {
242  vertices = vtxReco->vertices(cluster->tracks, bs); // attempt with config given reconstructor
243  }
244  TransientVertex singleFitVertex;
245  if (useVertexFitter) {
246  singleFitVertex = theAdaptiveFitter.vertex(cluster->tracks, cluster->seedPoint); //attempt with direct fitting
247  if (singleFitVertex.isValid())
248  vertices.push_back(singleFitVertex);
249  }
250 
251  // for each transient vertex state determine if a time can be measured and fill covariance
252  if (pv.covariance(3, 3) > 0.) {
253  for (auto &vtx : vertices) {
255  }
256  }
257 
258  for (std::vector<TransientVertex>::const_iterator v = vertices.begin(); v != vertices.end(); ++v) {
259  Measurement1D dlen = vdist.distance(pv, *v);
260  Measurement1D dlen2 = vdist2d.distance(pv, *v);
261 #ifdef VTXDEBUG
262  VTX vv(*v);
263  std::cout << "V chi2/n: " << v->normalisedChiSquared() << " ndof: " << v->degreesOfFreedom();
264  std::cout << " dlen: " << dlen.value() << " error: " << dlen.error() << " signif: " << dlen.significance();
265  std::cout << " dlen2: " << dlen2.value() << " error2: " << dlen2.error()
266  << " signif2: " << dlen2.significance();
267  std::cout << " pos: " << vv.position() << " error: " << vv.xError() << " " << vv.yError() << " " << vv.zError()
268  << std::endl;
269  std::cout << " time: " << vv.time() << " error: " << vv.tError() << std::endl;
270 #endif
272  std::vector<reco::TransientTrack> ts = v->originalTracks();
273  for (std::vector<reco::TransientTrack>::const_iterator i = ts.begin(); i != ts.end(); ++i) {
274  float w = v->trackWeight(*i);
275  if (w > 0.5)
276  dir += i->impactPointState().globalDirection();
277 #ifdef VTXDEBUG
278  std::cout << "\t[" << (*i).track().pt() << ": " << (*i).track().eta() << ", " << (*i).track().phi() << "], "
279  << w << std::endl;
280 #endif
281  }
282  GlobalPoint sv((*v).position().x(), (*v).position().y(), (*v).position().z());
283  float vscal = dir.unit().dot((sv - ppv).unit());
284  if (dlen.significance() > vertexMinDLenSig &&
285  ((vertexMinAngleCosine > 0) ? (vscal > vertexMinAngleCosine) : (vscal < vertexMinAngleCosine)) &&
286  v->normalisedChiSquared() < 10 && dlen2.significance() > vertexMinDLen2DSig) {
287  recoVertices->push_back(*v);
288 
289 #ifdef VTXDEBUG
290  std::cout << "ADDED" << std::endl;
291 #endif
292  }
293  }
294  }
295 #ifdef VTXDEBUG
296 
297  std::cout << "Final put " << recoVertices->size() << std::endl;
298 #endif
299  }
300 
301  event.put(std::move(recoVertices));
302 }
303 #endif
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:29
T getParameter(std::string const &) const
const double w
Definition: UKUtility.cc:23
void setBeamSpot(const reco::BeamSpot &beamSpot)
bool isValid() const
Make the ReferenceCountingProxy method to check validity public.
std::pair< std::vector< reco::TransientTrack >, GlobalPoint > nearTracks(const reco::TransientTrack &seed, const std::vector< reco::TransientTrack > &tracks, const reco::Vertex &primaryVertex) const
int numberOfValidHits() const
Definition: HitPattern.h:787
double error() const
Definition: Measurement1D.h:27
double timeExt() const
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
std::unique_ptr< VertexReconstructor > vtxReco
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:134
InputContainer::value_type TRK
constexpr bool isFinite(T x)
const Point & position() const
position
Definition: Vertex.h:113
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
std::unique_ptr< TracksClusteringFromDisplacedSeed > clusterizer
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:19
double pt() const
track transverse momentum
Definition: TrackBase.h:602
def pv(vc)
Definition: MetAnalyzer.py:7
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:596
Vector3DBase unit() const
Definition: Vector3DBase.h:54
double significance() const
Definition: Measurement1D.h:29
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:483
double dtErrorExt() const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double value() const
Definition: Measurement1D.h:25
void updateVertexTime(TransientVertex &vtx)
Definition: SVTimeHelpers.h:11
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:73
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
TemplatedInclusiveVertexFinder(const edm::ParameterSet &params)
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
reco::TransientTrack buildTT(edm::Handle< reco::TrackCollection > &tracks, edm::ESHandle< TransientTrackBuilder > &trackbuilder, unsigned int k)
Definition: TTHelpers.h:10
primaryVertex
hltOfflineBeamSpot for HLTMON
static void fillDescriptions(edm::ConfigurationDescriptions &cdesc)
bool isValid() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
double t() const
t coordinate
Definition: Vertex.h:121
Basic3DVector unit() const