CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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.size() != 0) {
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  FreeTrajectoryState primTraj = getTrajectory(primTrack);
90 
91  // get the secondary track with the max number of hits
92  unsigned short maxHits = 0; int indice=-1;
93  for( unsigned short i=0; i < secTracks.size(); i++) {
94  // Add references to daughters
95  unsigned short nhits = secTracks[i]->numberOfValidHits();
96  if( nhits > maxHits ) { maxHits = nhits; indice=i; }
97  }
98 
99  // Closest points
100  if(indice == -1) return false;
101 
102  ClosestApproachInRPhi* theApproach = closestApproach( primTrack, secTracks[indice]);
103  GlobalPoint crossing = theApproach->crossingPoint();
104  delete theApproach;
105 
106  // Create vertex (creation point)
107  // TODO: add error
109  crossing.y(),
110  crossing.z()),
111  reco::Vertex::Error(), 0.,0.,0);
112 
113  the_vertex.add(reco::TrackBaseRef( primTrack ), 1.0);
114  for( unsigned short i=0; i < secTracks.size(); i++) {
115  if( i==indice ) the_vertex.add(reco::TrackBaseRef( secTracks[i] ), 1.0);
116  else the_vertex.add(reco::TrackBaseRef( secTracks[i] ), 0.0);
117  }
118  return true;
119 }
120 
122  FreeTrajectoryState primTraj = getTrajectory(primTrack);
123  ClosestApproachInRPhi *theApproach = new ClosestApproachInRPhi();
124  FreeTrajectoryState secTraj = getTrajectory(secTrack);
125  bool status = theApproach->calculate(primTraj,secTraj);
126  if( status ) { return theApproach; }
127  else {
128  return NULL;
129  }
130 }
131 
132 bool NuclearVertexBuilder::isGoodSecondaryTrack( const reco::TrackRef& primTrack, const reco::TrackRef& secTrack ) const {
133  ClosestApproachInRPhi* theApproach = closestApproach(primTrack, secTrack);
134  bool result = false;
135  if(theApproach)
136  result = isGoodSecondaryTrack(secTrack, primTrack, theApproach->distance() , theApproach->crossingPoint());
137  delete theApproach;
138  return result;
139 }
140 
142  const reco::TrackRef& primTrack,
143  const double& distOfClosestApp,
144  const GlobalPoint& crossPoint ) const {
145  float TRACKER_RADIUS=129;
146  double pt2 = secTrack->pt();
147  double Dpt2 = secTrack->ptError();
148  double p1 = primTrack->p();
149  double Dp1 = primTrack->qoverpError()*p1*p1;
150  double p2 = secTrack->p();
151  double Dp2 = secTrack->qoverpError()*p2*p2;
152  std::cout << "1)" << distOfClosestApp << " < " << minDistFromPrim_ << " " << (distOfClosestApp < minDistFromPrim_) << "\n";
153  std::cout << "2)" << secTrack->normalizedChi2() << " < " << chi2Cut_ << " " << (secTrack->normalizedChi2() < chi2Cut_) << "\n";
154  std::cout << "3)" << crossPoint.perp() << " < " << TRACKER_RADIUS << " " << (crossPoint.perp() < TRACKER_RADIUS) << "\n";
155  std::cout << "4)" << (Dpt2/pt2) << " < " << DPtovPtCut_ << " " << ((Dpt2/pt2) < DPtovPtCut_) << "\n";
156  std::cout << "5)" << (p2-2*Dp2) << " < " << (p1+2*Dp1) << " " << ((p2-2*Dp2) < (p1+2*Dp1))<< "\n";
157  if( distOfClosestApp < minDistFromPrim_ &&
158  secTrack->normalizedChi2() < chi2Cut_ &&
159  crossPoint.perp() < TRACKER_RADIUS &&
160  (Dpt2/pt2) < DPtovPtCut_ &&
161  (p2-2*Dp2) < (p1+2*Dp1)) return true;
162  else return false;
163 }
164 
165 
167 
168  reco::TrackRef primTrack = (*(the_vertex.tracks_begin())).castTo<reco::TrackRef>();
169  ClosestApproachInRPhi *theApproach = closestApproach(primTrack, secTrack);
170  bool result = false;
171  if( theApproach ) {
172  GlobalPoint crp = theApproach->crossingPoint();
174  float dist = sqrt((crp.x()-vtx.x())*(crp.x()-vtx.x()) +
175  (crp.y()-vtx.y())*(crp.y()-vtx.y()) +
176  (crp.z()-vtx.z())*(crp.z()-vtx.z()));
177 
178  float distError = sqrt(the_vertex.xError()*the_vertex.xError() +
181 
182  //std::cout << "Distance between Additional track and vertex =" << dist << " +/- " << distError << std::endl;
183  // TODO : add condition on distance between last rechit of the primary and the first rec hit of the secondary
184 
185  result = (isGoodSecondaryTrack(secTrack, primTrack, theApproach->distance(), crp)
186  && dist-distError<minDistFromVtx_);
187  }
188  delete theApproach;
189  return result;
190 }
191 
193  std::vector<reco::TrackRef> allSecondary;
195  allSecondary.push_back( (*it).castTo<reco::TrackRef>() );
196  }
197  allSecondary.push_back( secTrack );
198  build( (*the_vertex.tracks_begin()).castTo<reco::TrackRef>(), allSecondary );
199 }
200 
202  std::vector<reco::TrackRef>& tC) const {
203 
204  // inspired from FinalTrackSelector (S. Wagner) modified by P. Janot
205  LogDebug("NuclearInteractionMaker") << "cleanTrackCollection number of input tracks : " << tC.size();
206  std::map<std::vector<reco::TrackRef>::const_iterator, std::vector<const TrackingRecHit*> > rh;
207 
208  // first remove bad quality tracks and create map
209  std::vector<bool> selected(tC.size(), false);
210  int i=0;
211  for (std::vector<reco::TrackRef>::const_iterator track=tC.begin(); track!=tC.end(); track++){
212  if( isGoodSecondaryTrack(primTrack, *track)) {
213  selected[i]=true;
214  trackingRecHit_iterator itB = (*track)->recHitsBegin();
215  trackingRecHit_iterator itE = (*track)->recHitsEnd();
216  for (trackingRecHit_iterator it = itB; it != itE; ++it) {
217  const TrackingRecHit* hit = &(**it);
218  rh[track].push_back(hit);
219  }
220  }
221  i++;
222  }
223 
224  // then remove duplicated tracks
225  i=-1;
226  for (std::vector<reco::TrackRef>::const_iterator track=tC.begin(); track!=tC.end(); track++){
227  i++;
228  int j=-1;
229  for (std::vector<reco::TrackRef>::const_iterator track2=tC.begin(); track2!=tC.end(); track2++){
230  j++;
231  if ((!selected[j])||(!selected[i]))continue;
232  if ((j<=i))continue;
233  int noverlap=0;
234  std::vector<const TrackingRecHit*>& iHits = rh[track];
235  for ( unsigned ih=0; ih<iHits.size(); ++ih ) {
236  const TrackingRecHit* it = iHits[ih];
237  if (it->isValid()){
238  std::vector<const TrackingRecHit*>& jHits = rh[track2];
239  for ( unsigned ih2=0; ih2<jHits.size(); ++ih2 ) {
240  const TrackingRecHit* jt = jHits[ih2];
241  if (jt->isValid()){
242  const TrackingRecHit* kt = jt;
243  if ( it->sharesInput(kt,TrackingRecHit::some) )noverlap++;
244  }
245  }
246  }
247  }
248  float fi=float(noverlap)/float((*track)->recHitsSize());
249  float fj=float(noverlap)/float((*track2)->recHitsSize());
250  if ((fi>shareFrac_)||(fj>shareFrac_)){
251  if (fi<fj){
252  selected[j]=false;
253  }else{
254  if (fi>fj){
255  selected[i]=false;
256  }else{
257  if ((*track)->normalizedChi2() > (*track2)->normalizedChi2()){selected[i]=false;}else{selected[j]=false;}
258  }//end fi > or = fj
259  }//end fi < fj
260  }//end got a duplicate
261  }//end track2 loop
262  }//end track loop
263 
264  std::vector< reco::TrackRef > newTrackColl;
265  i=0;
266  for (std::vector<reco::TrackRef>::const_iterator track=tC.begin(); track!=tC.end(); track++){
267  if( selected[i] ) newTrackColl.push_back( *track );
268  ++i;
269  }
270  tC = newTrackColl;
271 }
272 
274  std::vector<reco::TrackRef>& tC) const {
275  float totalEnergy=0;
276  for(size_t i=0; i< tC.size(); ++i) {
277  totalEnergy += tC[i]->p();
278  }
279  if( totalEnergy > primTrack->p()+0.1*primTrack->p() ) {
280  tC.pop_back();
281  checkEnergy(primTrack,tC);
282  }
283 }
#define LogDebug(id)
virtual char const * what() const
Definition: Exception.cc:97
void addSecondaryTrack(const reco::TrackRef &secTrack)
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:66
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
double zError() const
error on z
Definition: Vertex.h:105
Common base class.
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:61
reco::TransientTrack build(const reco::Track *p) const
T y() const
Definition: PV3DBase.h:57
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
void build(const reco::TrackRef &primaryTrack, std::vector< reco::TrackRef > &secondaryTrack)
#define NULL
Definition: scimark2.h:8
const Point & position() const
position
Definition: Vertex.h:93
const TransientTrackBuilder * theTransientTrackBuilder
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
virtual GlobalPoint crossingPoint() const
ClosestApproachInRPhi * closestApproach(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
T sqrt(T t)
Definition: SSEVec.h:28
virtual const char * what() const
T z() const
Definition: PV3DBase.h:58
tuple result
Definition: query.py:137
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
bool isCompatible(const reco::TrackRef &secTrack) const
int j
Definition: DBlmapReader.cc:9
void checkEnergy(const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
void FillVertexWithLastPrimHit(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
double p2[4]
Definition: TauolaWrapper.h:90
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
double xError() const
error on x
Definition: Vertex.h:101
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
bool isValid() const
double p1[4]
Definition: TauolaWrapper.h:89
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
bool FillVertexWithCrossingPoint(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
FreeTrajectoryState getTrajectory(const reco::TrackRef &track) const
tuple cout
Definition: gather_cfg.py:41
bool isGoodSecondaryTrack(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
const MagneticField * theMagField
tuple status
Definition: ntuplemaker.py:245
virtual float distance() const
T x() const
Definition: PV3DBase.h:56
void cleanTrackCollection(const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
double yError() const
error on y
Definition: Vertex.h:103
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:35
bool FillVertexWithAdaptVtxFitter(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)