CMS 3D CMS Logo

NuclearVertexBuilder.cc
Go to the documentation of this file.
7 
9 
10 void NuclearVertexBuilder::build( const reco::TrackRef& primTrack, std::vector<reco::TrackRef>& secTracks) {
11 
12  cleanTrackCollection(primTrack, secTracks);
13  std::sort(secTracks.begin(),secTracks.end(),cmpTracks());
14  checkEnergy(primTrack, secTracks);
15 
16  if( !secTracks.empty()) {
17  if( FillVertexWithAdaptVtxFitter(primTrack, secTracks) ) return;
18  else if( FillVertexWithCrossingPoint(primTrack, secTracks) ) return;
19  else FillVertexWithLastPrimHit( primTrack, secTracks);
20  }
21  else {
22  // if no secondary tracks : vertex position = position of last rechit of the primary track
23  FillVertexWithLastPrimHit( primTrack, secTracks);
24  }
25 }
26 
28 {
29  GlobalPoint position(track->vertex().x(),
30  track->vertex().y(),
31  track->vertex().z());
32 
33  GlobalVector momentum(track->momentum().x(),
34  track->momentum().y(),
35  track->momentum().z());
36 
38  track->charge(),theMagField);
39 
40  FreeTrajectoryState fts(gtp);
41 
42  return fts;
43 }
44 
45 void NuclearVertexBuilder::FillVertexWithLastPrimHit(const reco::TrackRef& primTrack, const std::vector<reco::TrackRef>& secTracks) {
46  LogDebug("NuclearInteractionMaker") << "Vertex build from the end point of the primary track";
47  the_vertex = reco::Vertex(reco::Vertex::Point(primTrack->outerPosition().x(),
48  primTrack->outerPosition().y(),
49  primTrack->outerPosition().z()),
50  reco::Vertex::Error(), 0.,0.,0);
51  the_vertex.add(reco::TrackBaseRef( primTrack ), 1.0);
52  for( unsigned short i=0; i != secTracks.size(); i++) {
53  the_vertex.add(reco::TrackBaseRef( secTracks[i] ), 0.0);
54  }
55 }
56 
57 bool NuclearVertexBuilder::FillVertexWithAdaptVtxFitter(const reco::TrackRef& primTrack, const std::vector<reco::TrackRef>& secTracks) {
58  std::vector<reco::TransientTrack> transientTracks;
59  transientTracks.push_back( theTransientTrackBuilder->build(primTrack));
60  // get the secondary track with the max number of hits
61  for( unsigned short i=0; i != secTracks.size(); i++ ) {
62  transientTracks.push_back( theTransientTrackBuilder->build( secTracks[i]) );
63  }
64  if( transientTracks.size() == 1 ) return false;
66  try {
67  TransientVertex tv = AVF.vertex(transientTracks);
69  }
70  catch(VertexException& exception){
71  // AdaptivevertexFitter does not work
72  LogDebug("NuclearInteractionMaker") << exception.what() << "\n";
73  return false;
74  }
75  catch( cms::Exception& exception){
76  // AdaptivevertexFitter does not work
77  LogDebug("NuclearInteractionMaker") << exception.what() << "\n";
78  return false;
79  }
80  if( the_vertex.isValid() ) {
81  LogDebug("NuclearInteractionMaker") << "Try to build from AdaptiveVertexFitter with " << the_vertex.tracksSize() << " tracks";
82  return true;
83  }
84  else return false;
85 }
86 
87 
88 bool NuclearVertexBuilder::FillVertexWithCrossingPoint(const reco::TrackRef& primTrack, const std::vector<reco::TrackRef>& secTracks) {
89  // get the secondary track with the max number of hits
90  unsigned short maxHits = 0; int indice=-1;
91  for( unsigned short i=0; i < secTracks.size(); i++) {
92  // Add references to daughters
93  unsigned short nhits = secTracks[i]->numberOfValidHits();
94  if( nhits > maxHits ) { maxHits = nhits; indice=i; }
95  }
96 
97  // Closest points
98  if(indice == -1) return false;
99 
100  ClosestApproachInRPhi* theApproach = closestApproach( primTrack, secTracks[indice]);
101  GlobalPoint crossing = theApproach->crossingPoint();
102  delete theApproach;
103 
104  // Create vertex (creation point)
105  // TODO: add error
107  crossing.y(),
108  crossing.z()),
109  reco::Vertex::Error(), 0.,0.,0);
110 
111  the_vertex.add(reco::TrackBaseRef( primTrack ), 1.0);
112  for( unsigned short i=0; i < secTracks.size(); i++) {
113  if( i==indice ) the_vertex.add(reco::TrackBaseRef( secTracks[i] ), 1.0);
114  else the_vertex.add(reco::TrackBaseRef( secTracks[i] ), 0.0);
115  }
116  return true;
117 }
118 
120  FreeTrajectoryState primTraj = getTrajectory(primTrack);
121  ClosestApproachInRPhi *theApproach = new ClosestApproachInRPhi();
122  FreeTrajectoryState secTraj = getTrajectory(secTrack);
123  bool status = theApproach->calculate(primTraj,secTraj);
124  if( status ) { return theApproach; }
125  else {
126  return nullptr;
127  }
128 }
129 
130 bool NuclearVertexBuilder::isGoodSecondaryTrack( const reco::TrackRef& primTrack, const reco::TrackRef& secTrack ) const {
131  ClosestApproachInRPhi* theApproach = closestApproach(primTrack, secTrack);
132  bool result = false;
133  if(theApproach)
134  result = isGoodSecondaryTrack(secTrack, primTrack, theApproach->distance() , theApproach->crossingPoint());
135  delete theApproach;
136  return result;
137 }
138 
140  const reco::TrackRef& primTrack,
141  const double& distOfClosestApp,
142  const GlobalPoint& crossPoint ) const {
143  float TRACKER_RADIUS=129;
144  double pt2 = secTrack->pt();
145  double Dpt2 = secTrack->ptError();
146  double p1 = primTrack->p();
147  double Dp1 = primTrack->qoverpError()*p1*p1;
148  double p2 = secTrack->p();
149  double Dp2 = secTrack->qoverpError()*p2*p2;
150  std::cout << "1)" << distOfClosestApp << " < " << minDistFromPrim_ << " " << (distOfClosestApp < minDistFromPrim_) << "\n";
151  std::cout << "2)" << secTrack->normalizedChi2() << " < " << chi2Cut_ << " " << (secTrack->normalizedChi2() < chi2Cut_) << "\n";
152  std::cout << "3)" << crossPoint.perp() << " < " << TRACKER_RADIUS << " " << (crossPoint.perp() < TRACKER_RADIUS) << "\n";
153  std::cout << "4)" << (Dpt2/pt2) << " < " << DPtovPtCut_ << " " << ((Dpt2/pt2) < DPtovPtCut_) << "\n";
154  std::cout << "5)" << (p2-2*Dp2) << " < " << (p1+2*Dp1) << " " << ((p2-2*Dp2) < (p1+2*Dp1))<< "\n";
155  if( distOfClosestApp < minDistFromPrim_ &&
156  secTrack->normalizedChi2() < chi2Cut_ &&
157  crossPoint.perp() < TRACKER_RADIUS &&
158  (Dpt2/pt2) < DPtovPtCut_ &&
159  (p2-2*Dp2) < (p1+2*Dp1)) return true;
160  else return false;
161 }
162 
163 
165 
166  reco::TrackRef primTrack = (*(the_vertex.tracks_begin())).castTo<reco::TrackRef>();
167  ClosestApproachInRPhi *theApproach = closestApproach(primTrack, secTrack);
168  bool result = false;
169  if( theApproach ) {
170  GlobalPoint crp = theApproach->crossingPoint();
172  float dist = sqrt((crp.x()-vtx.x())*(crp.x()-vtx.x()) +
173  (crp.y()-vtx.y())*(crp.y()-vtx.y()) +
174  (crp.z()-vtx.z())*(crp.z()-vtx.z()));
175 
176  float distError = sqrt(the_vertex.xError()*the_vertex.xError() +
179 
180  //std::cout << "Distance between Additional track and vertex =" << dist << " +/- " << distError << std::endl;
181  // TODO : add condition on distance between last rechit of the primary and the first rec hit of the secondary
182 
183  result = (isGoodSecondaryTrack(secTrack, primTrack, theApproach->distance(), crp)
184  && dist-distError<minDistFromVtx_);
185  }
186  delete theApproach;
187  return result;
188 }
189 
191  std::vector<reco::TrackRef> allSecondary;
193  allSecondary.push_back( (*it).castTo<reco::TrackRef>() );
194  }
195  allSecondary.push_back( secTrack );
196  build( (*the_vertex.tracks_begin()).castTo<reco::TrackRef>(), allSecondary );
197 }
198 
200  std::vector<reco::TrackRef>& tC) const {
201 
202  // inspired from FinalTrackSelector (S. Wagner) modified by P. Janot
203  LogDebug("NuclearInteractionMaker") << "cleanTrackCollection number of input tracks : " << tC.size();
204  std::map<std::vector<reco::TrackRef>::const_iterator, std::vector<const TrackingRecHit*> > rh;
205 
206  // first remove bad quality tracks and create map
207  std::vector<bool> selected(tC.size(), false);
208  int i=0;
209  for (std::vector<reco::TrackRef>::const_iterator track=tC.begin(); track!=tC.end(); track++){
210  if( isGoodSecondaryTrack(primTrack, *track)) {
211  selected[i]=true;
212  trackingRecHit_iterator itB = (*track)->recHitsBegin();
213  trackingRecHit_iterator itE = (*track)->recHitsEnd();
214  for (trackingRecHit_iterator it = itB; it != itE; ++it) {
215  const TrackingRecHit* hit = &(**it);
216  rh[track].push_back(hit);
217  }
218  }
219  i++;
220  }
221 
222  // then remove duplicated tracks
223  i=-1;
224  for (std::vector<reco::TrackRef>::const_iterator track=tC.begin(); track!=tC.end(); track++){
225  i++;
226  int j=-1;
227  for (std::vector<reco::TrackRef>::const_iterator track2=tC.begin(); track2!=tC.end(); track2++){
228  j++;
229  if ((!selected[j])||(!selected[i]))continue;
230  if ((j<=i))continue;
231  int noverlap=0;
232  std::vector<const TrackingRecHit*>& iHits = rh[track];
233  for ( unsigned ih=0; ih<iHits.size(); ++ih ) {
234  const TrackingRecHit* it = iHits[ih];
235  if (it->isValid()){
236  std::vector<const TrackingRecHit*>& jHits = rh[track2];
237  for ( unsigned ih2=0; ih2<jHits.size(); ++ih2 ) {
238  const TrackingRecHit* jt = jHits[ih2];
239  if (jt->isValid()){
240  const TrackingRecHit* kt = jt;
241  if ( it->sharesInput(kt,TrackingRecHit::some) )noverlap++;
242  }
243  }
244  }
245  }
246  float fi=float(noverlap)/float((*track)->recHitsSize());
247  float fj=float(noverlap)/float((*track2)->recHitsSize());
248  if ((fi>shareFrac_)||(fj>shareFrac_)){
249  if (fi<fj){
250  selected[j]=false;
251  }else{
252  if (fi>fj){
253  selected[i]=false;
254  }else{
255  if ((*track)->normalizedChi2() > (*track2)->normalizedChi2()){selected[i]=false;}else{selected[j]=false;}
256  }//end fi > or = fj
257  }//end fi < fj
258  }//end got a duplicate
259  }//end track2 loop
260  }//end track loop
261 
262  std::vector< reco::TrackRef > newTrackColl;
263  i=0;
264  for (std::vector<reco::TrackRef>::const_iterator track=tC.begin(); track!=tC.end(); track++){
265  if( selected[i] ) newTrackColl.push_back( *track );
266  ++i;
267  }
268  tC = newTrackColl;
269 }
270 
272  std::vector<reco::TrackRef>& tC) const {
273  float totalEnergy=0;
274  for(size_t i=0; i< tC.size(); ++i) {
275  totalEnergy += tC[i]->p();
276  }
277  if( totalEnergy > primTrack->p()+0.1*primTrack->p() ) {
278  tC.pop_back();
279  checkEnergy(primTrack,tC);
280  }
281 }
#define LogDebug(id)
void addSecondaryTrack(const reco::TrackRef &secTrack)
float distance() const override
T perp() const
Definition: PV3DBase.h:72
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
double zError() const
error on z
Definition: Vertex.h:123
Common base class.
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:68
reco::TransientTrack build(const reco::Track *p) const
T y() const
Definition: PV3DBase.h:63
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
void build(const reco::TrackRef &primaryTrack, std::vector< reco::TrackRef > &secondaryTrack)
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
char const * what() const override
Definition: Exception.cc:141
const Point & position() const
position
Definition: Vertex.h:109
const TransientTrackBuilder * theTransientTrackBuilder
ClosestApproachInRPhi * closestApproach(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
bool isCompatible(const reco::TrackRef &secTrack) const
void checkEnergy(const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
void FillVertexWithLastPrimHit(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
double p2[4]
Definition: TauolaWrapper.h:90
const char * what() const override
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
double xError() const
error on x
Definition: Vertex.h:119
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool isValid() const
double p1[4]
Definition: TauolaWrapper.h:89
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
bool FillVertexWithCrossingPoint(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
FreeTrajectoryState getTrajectory(const reco::TrackRef &track) const
bool isGoodSecondaryTrack(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
const MagneticField * theMagField
GlobalPoint crossingPoint() const override
T x() const
Definition: PV3DBase.h:62
void cleanTrackCollection(const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
double yError() const
error on y
Definition: Vertex.h:121
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:71
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
bool FillVertexWithAdaptVtxFitter(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)