CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFDisplacedVertexCandidateFinder.cc
Go to the documentation of this file.
2 
5 
7 
8 using namespace std;
9 using namespace reco;
10 
11 //for debug only
12 //#define PFLOW_DEBUG
13 
15  vertexCandidates_( new PFDisplacedVertexCandidateCollection ),
16  dcaCut_(1000),
17  primaryVertexCut2_(0.0),
18  dcaPInnerHitCut2_(1000.0),
19  vertexCandidatesSize_(50),
20  debug_(false) {
21 
23  theMinimum_ = theMinimum;
24 }
25 
26 
28 
29 #ifdef PFLOW_DEBUG
30  if(debug_)
31  cout<<"~PFDisplacedVertexCandidateFinder - number of remaining elements: "
32  <<eventTracks_.size()<<endl;
33 #endif
34 
35 }
36 
37 // Set the imput collection of tracks and calculate their
38 // trajectory parameters the Global Trajectory Parameters
39 void
41  const MagneticField* magField) {
42 
43  magField_ = magField;
44  trackMask_.clear();
45  trackMask_.resize(trackh->size());
47  eventTrackTrajectories_.resize(trackh->size());
48 
49  for(unsigned i=0;i<trackMask_.size(); i++) trackMask_[i] = true;
50 
51  eventTracks_.clear();
52  if(trackh.isValid()) {
53  for(unsigned i=0;i<trackh->size(); i++) {
54 
55  TrackRef tref( trackh, i);
56  TrackBaseRef tbref(tref);
57 
58  if( !isSelected( tbref ) ) {
59  trackMask_[i] = false;
60  continue;
61  }
62 
63  const Track* trk = tref.get();
64  eventTracks_.push_back( tbref );
66  }
67  }
68 }
69 
70 
71 // -------- Main function which find vertices -------- //
72 
73 void
75 
76  if (debug_) cout << "========= Start Finding Displaced Vertex Candidates =========" << endl;
77 
78  // The vertexCandidates have not been passed to the event, and need to be cleared
79  if(vertexCandidates_.get() )vertexCandidates_->clear();
80  else
82 
84  for(IE ie = eventTracks_.begin();
85  ie != eventTracks_.end();) {
86 
87  // Run the recursive procedure to find all tracks link together
88  // In one blob called Candidate
89 
90  PFDisplacedVertexCandidate tempVertexCandidate;
91 
92  ie = associate( eventTracks_.end() , ie, tempVertexCandidate);
93 
94  // Build remaining links in current block
95  if(tempVertexCandidate.isValid()) {
96  packLinks( tempVertexCandidate );
97  vertexCandidates_->push_back( tempVertexCandidate );
98  }
99  }
100 
101  if(debug_) cout << "========= End Finding Displaced Vertex Candidates =========" << endl;
102 
103 }
104 
105 
106 
107 
110  PFDisplacedVertexCandidate& tempVertexCandidate) {
111 
112 
113 #ifdef PFLOW_DEBUG
114  if(debug_ ) cout<<"== Start the association procedure =="<<endl;
115 #endif
116 
117  if( last!= eventTracks_.end() ) {
118  double dist = -1;
119  GlobalPoint P(0,0,0);
121  link( (*last), (*next), dist, P, linktest);
122 
123 
124  if(dist<-0.5) {
125 #ifdef PFLOW_DEBUG
126  if(debug_ ) cout<<"link failed"<<endl;
127 #endif
128  return ++next; // association failed
129  }
130  else {
131  // add next element to the current pflowblock
132  tempVertexCandidate.addElement( (*next));
133  trackMask_[(*next).key()] = false;
134 #ifdef PFLOW_DEBUG
135  if(debug_ )
136  cout<<"link parameters "
137  << " *next = " << (*next).key()
138  << " *last = " << (*last).key()
139  << " dist = " << dist
140  << " P.x = " << P.x()
141  << " P.y = " << P.y()
142  << " P.z = " << P.z()
143  << endl;
144 #endif
145  }
146  }
147  else {
148  // add next element to this eflowblock
149 #ifdef PFLOW_DEBUG
150  if(debug_ ) cout<<"adding to block element "
151  << (*next).key()
152  <<endl;
153 #endif
154  tempVertexCandidate.addElement( (*next));
155  trackMask_[(*next).key()] = false;
156 
157 
158  }
159 
160  // recursive call: associate next and other unused elements
161 #ifdef PFLOW_DEBUG
162  if(debug_ ) {
163  for(unsigned i=0; i<trackMask_.size(); i++)
164  cout << " Mask[" << i << "] = " << trackMask_[i];
165  cout << "" << endl;
166  }
167 #endif
168 
169  for(IE ie = eventTracks_.begin();
170  ie != eventTracks_.end();) {
171 
172  if( ie == last || ie == next ) { ++ie; continue;}
173 
174  // *ie already included to a block
175  if( !trackMask_[(*ie).key()] ) { ++ie; continue;}
176 
177 #ifdef PFLOW_DEBUG
178  if(debug_ ) cout<<"calling associate "
179  << (*next).key()
180  <<" & "
181  << (*ie).key()
182  <<endl;
183 #endif
184  ie = associate(next, ie, tempVertexCandidate);
185 
186  }
187 
188 #ifdef PFLOW_DEBUG
189  if(debug_ ) {
190  cout<<"**** removing element "<<endl;
191  }
192 #endif
193 
194  IE iteratorToNextFreeElement = eventTracks_.erase( next );
195 
196 #ifdef PFLOW_DEBUG
197  if(debug_ ) cout<< "== End the association procedure ==" <<endl;
198 #endif
199 
200  return iteratorToNextFreeElement;
201 }
202 
203 
204 
205 void
207  const TrackBaseRef& el2,
208  double& dist, GlobalPoint& P,
210  ){
211 
212  if ( fabs(el1->eta()-el2->eta()) > 1) {dist = -1; return;}
213  if ( el1->pt()>2 && el2->pt()>2 && fabs(el1->phi()-el2->phi()) > 1) {dist = -1; return;}
214 
217 
218 
219  // Closest approach
220  theMinimum_.calculate(gt1,gt2);
221 
222  // Fill the parameters
223  dist = theMinimum_.distance();
225 
226  vertexLinkTest = PFDisplacedVertexCandidate::LINKTEST_DCA; //rechit by default
227 
228 
229  // Check if the link test fails
230  if (dist > dcaCut_) {dist = -1; return;}
231 
232  // Check if the closses approach point is too close to the primary vertex/beam pipe
233  double rho2 = P.x()*P.x()+P.y()*P.y();
234 
235  if ( rho2 < primaryVertexCut2_) {dist = -1; return;}
236 
237  // Check if the closses approach point is out of range of the tracker
238 
239  double tob_rho_limit2 = 10000;
240  double tec_z_limit = 270;
241 
242  if ( rho2 > tob_rho_limit2) {dist = -1; return;}
243  if ( fabs(P.z()) > tec_z_limit) {dist = -1; return;}
244 
245  return;
246 
247 
248  /*
249  // Check if the inner hit of one of the tracks is not too far away from the vertex
250  double dx = el1->innerPosition().x()-P.x();
251  double dy = el1->innerPosition().y()-P.y();
252  double dz = el1->innerPosition().z()-P.z();
253  double dist2 = dx*dx+dy*dy+dz*dz;
254 
255  if (dist2 > dcaPInnerHitCut2_) {
256  dist = -1;
257  if (debug_) cout << "track " << el1.key() << " dist to vertex " << sqrt(dist2) << endl;
258  return;
259  }
260 
261  dx = el2->innerPosition().x()-P.x();
262  dy = el2->innerPosition().y()-P.y();
263  dz = el2->innerPosition().z()-P.z();
264  dist2 = dx*dx+dy*dy+dz*dz;
265 
266  if (dist2 > dcaPInnerHitCut2_) {
267  dist = -1;
268  if (debug_) cout << "track " << el2.key() << " dist to vertex " << sqrt(dist2) << endl;
269  return;
270  }
271  */
272 
273 }
274 
275 
276 // Build up a matrix of all the links between different tracks in
277 // In the Candidate
278 void
280 
281 
282  const vector < TrackBaseRef >& els = vertexCandidate.elements();
283 
284  //First Loop: update all link data
285  for( unsigned i1=0; i1<els.size(); i1++ ) {
286  for( unsigned i2=i1+1; i2<els.size(); i2++ ) {
287 
288  // no reflexive link
289  if( i1==i2 ) continue;
290 
291  double dist = -1;
292  GlobalPoint P(0,0,0);
294 
295  link( els[i1], els[i2], dist, P, linktest);
296 
297 
298 #ifdef PFLOW_DEBUG
299  if( debug_ )
300  cout << "Setting link between elements " << i1 << " key " << els[i1].key()
301  << " and " << i2 << " key " << els[i2].key()
302  << " of dist =" << dist << " computed from link test "
303  << linktest << endl;
304 #endif
305 
306  if(dist >-0.5) vertexCandidate.setLink( i1, i2, dist, P, linktest );
307  }
308  }
309 
310 }
311 
312 
313 
314 // --------------- TOOLS -------------- //
315 
316 
317 
318 // This tool is a copy from VZeroFinder
321 (const Track* track) const
322 {
323 
324  const GlobalPoint position(track->vx(),
325  track->vy(),
326  track->vz());
327 
328  const GlobalVector momentum(track->momentum().x(),
329  track->momentum().y(),
330  track->momentum().z());
331 
333  track->charge(),magField_);
334 
335  return gtp;
336 }
337 
338 
339 // This tool allow to pre-select the track for the displaced vertex finder
340 bool
342 
343  double nChi2 = trackref->normalizedChi2();
344  double pt = trackref->pt();
345  double dpt = trackref->ptError();
346 
347  double pt_error = dpt/pt*100;
348 
349  if (debug_) cout << " PFDisplacedVertexFinder: PFrecTrack->Track Pt= "
350  << pt << " dPt/Pt = " << pt_error << "% nChi2 = " << nChi2 << endl;
351  if (nChi2 > nChi2_max_ || pt < pt_min_){
352 
353  if (debug_) cout << " PFBlockAlgo: skip badly measured or low pt track"
354  << " nChi2_cut = " << 5
355  << " pt_cut = " << 0.2 << endl;
356  return false;
357  }
358 
359  return true;
360 }
361 
362 
363 
364 
365 
366 
367 
368 
369 
370 
371 
372 ostream& operator<<(std::ostream& out, const PFDisplacedVertexCandidateFinder& a) {
373  if(! out) return out;
374 
375  out<<"====== Particle Flow Block Algorithm ======= ";
376  out<<endl;
377  out<<"number of unassociated elements : "<<a.eventTracks_.size()<<endl;
378  out<<endl;
379 
381  ie != a.eventTracks_.end(); ie++) {
382 
383  double pt = (*ie).get()->pt();
384 
385  math::XYZPoint Pi = (*ie).get()->innerPosition();
386  math::XYZPoint Po = (*ie).get()->outerPosition();
387 
388  double innermost_radius = sqrt(Pi.x()*Pi.x() + Pi.y()*Pi.y() + Pi.z()*Pi.z());
389  double outermost_radius = sqrt(Po.x()*Po.x() + Po.y()*Po.y() + Po.z()*Po.z());
390  double innermost_rho = sqrt(Pi.x()*Pi.x() + Pi.y()*Pi.y());
391  double outermost_rho = sqrt(Po.x()*Po.x() + Po.y()*Po.y());
392 
393  out<<"ie = " << (*ie).key()
394  <<" pt = " << pt
395  <<" innermost hit radius = " << innermost_radius << " rho = " << innermost_rho
396  <<" outermost hit radius = " << outermost_radius << " rho = " << outermost_rho
397  <<endl;
398  }
399 
400 
401  const std::auto_ptr< reco::PFDisplacedVertexCandidateCollection >& vertexCandidates
402  = a.vertexCandidates();
403 
404  if(!vertexCandidates.get() ) {
405  out<<"vertexCandidates already transfered"<<endl;
406  }
407  else {
408  out<<"number of vertexCandidates : "<<vertexCandidates->size()<<endl;
409  out<<endl;
410 
411 
412  for(PFDisplacedVertexCandidateFinder::IBC ib=vertexCandidates->begin();
413  ib != vertexCandidates->end(); ib++)
414  ib->Dump();
415 
416 
417 
418  }
419 
420  return out;
421 }
virtual float distance() const
void setInput(const edm::Handle< reco::TrackCollection > &trackh, const MagneticField *magField)
const double Pi
void addElement(const TrackBaseRef)
add a track Reference to the current Candidate
void link(const reco::TrackBaseRef &el1, const reco::TrackBaseRef &el2, double &dist, GlobalPoint &P, reco::PFDisplacedVertexCandidate::VertexLinkTest &linktest)
Check whether 2 elements are linked and fill the link parameters.
int i
Definition: DBlmapReader.cc:9
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
std::vector< PFDisplacedVertexCandidate > PFDisplacedVertexCandidateCollection
collection of PFDisplacedVertexCandidate objects
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:150
A block of tracks linked together.
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:111
GlobalTrajectoryParameters getGlobalTrajectoryParameters(const reco::Track *) const
--—— TOOLS --—— //
virtual GlobalPoint crossingPoint() const
IE associate(IE next, IE last, reco::PFDisplacedVertexCandidate &tempVertexCandidate)
--—— Different steps of the finder algorithm --—— ///
T y() const
Definition: PV3DBase.h:57
std::vector< GlobalTrajectoryParameters > eventTrackTrajectories_
#define P
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:139
tuple els
Definition: asciidump.py:420
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:187
std::auto_ptr< reco::PFDisplacedVertexCandidateCollection > vertexCandidates_
-----— Members -----— ///
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
std::list< reco::TrackBaseRef >::const_iterator IEC
void setLink(unsigned i1, unsigned i2, const float dist, const GlobalPoint &dcaPoint, const VertexLinkTest test=LINKTEST_DCA)
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
T sqrt(T t)
Definition: SSEVec.h:28
double dcaCut_
–— Algo parameters for the vertex finder -— ///
double pt() const
track transverse momentum
Definition: TrackBase.h:131
T z() const
Definition: PV3DBase.h:58
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:194
bool isValid() const
A Vertex Candidate is valid if it has at least two tracks.
bool isValid() const
Definition: HandleBase.h:76
size_t key() const
Definition: RefToBase.h:228
tuple out
Definition: dbtoconf.py:99
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
reco::PFDisplacedVertexCandidateCollection::const_iterator IBC
const std::auto_ptr< reco::PFDisplacedVertexCandidateCollection > & vertexCandidates() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
bool goodPtResolution(const reco::TrackBaseRef &trackref) const
Quality Criterion on the Pt resolution to select a Track.
bool isSelected(const reco::TrackBaseRef &trackref)
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
double a
Definition: hdecay.h:121
void packLinks(reco::PFDisplacedVertexCandidate &vertexCandidate)
void findDisplacedVertexCandidates()
-----— Main function which find vertices -----— ///
double primaryVertexCut2_
Do not reconstruct vertices wich are too close to the beam pipe.
std::list< reco::TrackBaseRef > eventTracks_
The track refs.
tuple cout
Definition: gather_cfg.py:41
int charge() const
track electric charge
Definition: TrackBase.h:113
const std::vector< TrackBaseRef > & elements() const
std::list< reco::TrackBaseRef >::iterator IE
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:241
T x() const
Definition: PV3DBase.h:56
bool debug_
if true, debug printouts activated
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143