CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
NuclearVertexBuilder Class Reference

#include <NuclearVertexBuilder.h>

Classes

class  cmpTracks
 

Public Member Functions

void addSecondaryTrack (const reco::TrackRef &secTrack)
 
void build (const reco::TrackRef &primaryTrack, std::vector< reco::TrackRef > &secondaryTrack)
 
ClosestApproachInRPhiclosestApproach (const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
 
reco::Vertex getVertex () const
 
bool isCompatible (const reco::TrackRef &secTrack) const
 
 NuclearVertexBuilder (const MagneticField *mag, const TransientTrackBuilder *transientTkBuilder, const edm::ParameterSet &iConfig)
 

Private Member Functions

void checkEnergy (const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
 
void cleanTrackCollection (const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
 
bool FillVertexWithAdaptVtxFitter (const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
 
bool FillVertexWithCrossingPoint (const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
 
void FillVertexWithLastPrimHit (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
 
bool isGoodSecondaryTrack (const reco::TrackRef &secTrack, const reco::TrackRef &primTrack, const double &distOfClosestApp, const GlobalPoint &crossPoint) const
 

Private Attributes

double chi2Cut_
 
double DPtovPtCut_
 
double minDistFromPrim_
 
double minDistFromVtx_
 
double shareFrac_
 
reco::Vertex the_vertex
 
const MagneticFieldtheMagField
 
const TransientTrackBuildertheTransientTrackBuilder
 

Detailed Description

Definition at line 16 of file NuclearVertexBuilder.h.

Constructor & Destructor Documentation

NuclearVertexBuilder::NuclearVertexBuilder ( const MagneticField mag,
const TransientTrackBuilder transientTkBuilder,
const edm::ParameterSet iConfig 
)
inline

Definition at line 19 of file NuclearVertexBuilder.h.

References build().

19  :
20  theMagField(mag),
21  theTransientTrackBuilder(transientTkBuilder),
22  minDistFromPrim_( iConfig.getParameter<double>("minDistFromPrimary") ),
23  chi2Cut_(iConfig.getParameter<double>("chi2Cut")),
24  DPtovPtCut_(iConfig.getParameter<double>("DPtovPtCut")),
25  minDistFromVtx_(iConfig.getParameter<double>("minDistFromVtx")),
26  shareFrac_(iConfig.getParameter<double>("shareFrac")){}
T getParameter(std::string const &) const
const TransientTrackBuilder * theTransientTrackBuilder
const MagneticField * theMagField

Member Function Documentation

void NuclearVertexBuilder::addSecondaryTrack ( const reco::TrackRef secTrack)

Definition at line 190 of file NuclearVertexBuilder.cc.

References build(), the_vertex, reco::Vertex::tracks_begin(), and reco::Vertex::tracks_end().

Referenced by getVertex().

190  {
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 }
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
void build(const reco::TrackRef &primaryTrack, std::vector< reco::TrackRef > &secondaryTrack)
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
void NuclearVertexBuilder::build ( const reco::TrackRef primaryTrack,
std::vector< reco::TrackRef > &  secondaryTrack 
)

Definition at line 10 of file NuclearVertexBuilder.cc.

References checkEnergy(), cleanTrackCollection(), FillVertexWithAdaptVtxFitter(), FillVertexWithCrossingPoint(), and FillVertexWithLastPrimHit().

Referenced by addSecondaryTrack(), and NuclearVertexBuilder().

10  {
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 }
void checkEnergy(const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
void FillVertexWithLastPrimHit(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
bool FillVertexWithCrossingPoint(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
void cleanTrackCollection(const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
bool FillVertexWithAdaptVtxFitter(const reco::TrackRef &primTrack, const std::vector< reco::TrackRef > &secTracks)
void NuclearVertexBuilder::checkEnergy ( const reco::TrackRef primTrack,
std::vector< reco::TrackRef > &  tC 
) const
private

Definition at line 266 of file NuclearVertexBuilder.cc.

References mps_fire::i.

Referenced by build(), and getVertex().

267  {
268  float totalEnergy=0;
269  for(size_t i=0; i< tC.size(); ++i) {
270  totalEnergy += tC[i]->p();
271  }
272  if( totalEnergy > primTrack->p()+0.1*primTrack->p() ) {
273  tC.pop_back();
274  checkEnergy(primTrack,tC);
275  }
276 }
void checkEnergy(const reco::TrackRef &primTrack, std::vector< reco::TrackRef > &tC) const
void NuclearVertexBuilder::cleanTrackCollection ( const reco::TrackRef primTrack,
std::vector< reco::TrackRef > &  tC 
) const
private

Definition at line 199 of file NuclearVertexBuilder.cc.

References objects.autophobj::float, mps_fire::i, isGoodSecondaryTrack(), TrackingRecHit::isValid(), LogDebug, shareFrac_, TrackingRecHit::sharesInput(), TrackingRecHit::some, and HiIsolationCommonParameters_cff::track.

Referenced by build(), and getVertex().

200  {
201 
202  // inspired from FinalTrackSelector (S. Wagner) modified by P. Janot
203  LogDebug("NuclearInteractionMaker") << "cleanTrackCollection number of input tracks : " << tC.size();
204  std::map<reco::TrackRef const*, 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(auto const& track : tC) {
210  if( isGoodSecondaryTrack(primTrack, track)) {
211  selected[i]=true;
212  for(auto const& hit : track->recHits()) rh[&track].push_back(hit);
213  }
214  i++;
215  }
216 
217  // then remove duplicated tracks
218  i=-1;
219  for(auto const& track : tC) {
220  i++;
221  int j=-1;
222  for(auto const& track2 : tC) {
223  j++;
224  if ((!selected[j])||(!selected[i]))continue;
225  if ((j<=i))continue;
226  int noverlap=0;
227  std::vector<const TrackingRecHit*>& iHits = rh[&track];
228  for ( unsigned ih=0; ih<iHits.size(); ++ih ) {
229  const TrackingRecHit* it = iHits[ih];
230  if (it->isValid()){
231  std::vector<const TrackingRecHit*>& jHits = rh[&track2];
232  for ( unsigned ih2=0; ih2<jHits.size(); ++ih2 ) {
233  const TrackingRecHit* jt = jHits[ih2];
234  if (jt->isValid()){
235  const TrackingRecHit* kt = jt;
236  if ( it->sharesInput(kt,TrackingRecHit::some) )noverlap++;
237  }
238  }
239  }
240  }
241  float fi=float(noverlap)/float(track->recHitsSize());
242  float fj=float(noverlap)/float(track2->recHitsSize());
243  if ((fi>shareFrac_)||(fj>shareFrac_)){
244  if (fi<fj){
245  selected[j]=false;
246  }else{
247  if (fi>fj){
248  selected[i]=false;
249  }else{
250  if (track->normalizedChi2() > track2->normalizedChi2()){selected[i]=false;}else{selected[j]=false;}
251  }//end fi > or = fj
252  }//end fi < fj
253  }//end got a duplicate
254  }//end track2 loop
255  }//end track loop
256 
257  std::vector< reco::TrackRef > newTrackColl;
258  i=0;
259  for(auto const& track : tC) {
260  if( selected[i] ) newTrackColl.push_back(track);
261  ++i;
262  }
263  tC = newTrackColl;
264 }
#define LogDebug(id)
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
bool isValid() const
bool isGoodSecondaryTrack(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
ClosestApproachInRPhi * NuclearVertexBuilder::closestApproach ( const reco::TrackRef primTrack,
const reco::TrackRef secTrack 
) const

Definition at line 119 of file NuclearVertexBuilder.cc.

References ClosestApproachInRPhi::calculate(), getTrajectory(), and mps_update::status.

Referenced by FillVertexWithCrossingPoint(), NuclearInteractionEDProducer::findAdditionalSecondaryTracks(), getVertex(), isCompatible(), and isGoodSecondaryTrack().

119  {
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 }
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
FreeTrajectoryState getTrajectory(const reco::TrackRef &track) const
bool NuclearVertexBuilder::FillVertexWithAdaptVtxFitter ( const reco::TrackRef primTrack,
const std::vector< reco::TrackRef > &  secTracks 
)
private

Definition at line 57 of file NuclearVertexBuilder.cc.

References TransientTrackBuilder::build(), cppFunctionSkipper::exception, mps_fire::i, reco::Vertex::isValid(), LogDebug, the_vertex, theTransientTrackBuilder, reco::Vertex::tracksSize(), AdaptiveVertexFitter::vertex(), VertexException::what(), and cms::Exception::what().

Referenced by build(), and getVertex().

57  {
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 }
#define LogDebug(id)
Common base class.
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:68
reco::TransientTrack build(const reco::Track *p) const
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
char const * what() const override
Definition: Exception.cc:103
const TransientTrackBuilder * theTransientTrackBuilder
const char * what() const override
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:71
bool NuclearVertexBuilder::FillVertexWithCrossingPoint ( const reco::TrackRef primTrack,
const std::vector< reco::TrackRef > &  secTracks 
)
private

Definition at line 88 of file NuclearVertexBuilder.cc.

References reco::Vertex::add(), closestApproach(), ClosestApproachInRPhi::crossingPoint(), mps_fire::i, NuclearSeed_cfi::maxHits, nhits, the_vertex, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by build(), and getVertex().

88  {
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 }
T y() const
Definition: PV3DBase.h:63
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
ClosestApproachInRPhi * closestApproach(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
T z() const
Definition: PV3DBase.h:64
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
GlobalPoint crossingPoint() const override
T x() const
Definition: PV3DBase.h:62
void NuclearVertexBuilder::FillVertexWithLastPrimHit ( const reco::TrackRef primTrack,
const std::vector< reco::TrackRef > &  secTracks 
)
private

Definition at line 45 of file NuclearVertexBuilder.cc.

References reco::Vertex::add(), mps_fire::i, LogDebug, and the_vertex.

Referenced by build(), and getVertex().

45  {
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 }
#define LogDebug(id)
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
FreeTrajectoryState NuclearVertexBuilder::getTrajectory ( const reco::TrackRef track) const
private

Definition at line 27 of file NuclearVertexBuilder.cc.

References position, and theMagField.

Referenced by closestApproach(), and getVertex().

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 }
static int position[264][3]
Definition: ReadPGInfo.cc:509
const MagneticField * theMagField
reco::Vertex NuclearVertexBuilder::getVertex ( ) const
inline
bool NuclearVertexBuilder::isCompatible ( const reco::TrackRef secTrack) const

Definition at line 164 of file NuclearVertexBuilder.cc.

References closestApproach(), isGoodSecondaryTrack(), minDistFromVtx_, reco::Vertex::position(), mps_fire::result, mathSSE::sqrt(), the_vertex, reco::Vertex::tracks_begin(), extraflags_cff::vtx, PV3DBase< T, PVType, FrameType >::x(), reco::Vertex::xError(), PV3DBase< T, PVType, FrameType >::y(), reco::Vertex::yError(), PV3DBase< T, PVType, FrameType >::z(), and reco::Vertex::zError().

Referenced by getVertex().

164  {
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 }
double zError() const
error on z
Definition: Vertex.h:123
T y() const
Definition: PV3DBase.h:63
const Point & position() const
position
Definition: Vertex.h:109
ClosestApproachInRPhi * closestApproach(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
double xError() const
error on x
Definition: Vertex.h:119
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool isGoodSecondaryTrack(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
T x() const
Definition: PV3DBase.h:62
double yError() const
error on y
Definition: Vertex.h:121
bool NuclearVertexBuilder::isGoodSecondaryTrack ( const reco::TrackRef primTrack,
const reco::TrackRef secTrack 
) const
private

Definition at line 130 of file NuclearVertexBuilder.cc.

References closestApproach(), ClosestApproachInRPhi::crossingPoint(), ClosestApproachInRPhi::distance(), and mps_fire::result.

Referenced by cleanTrackCollection(), getVertex(), and isCompatible().

130  {
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 }
float distance() const override
ClosestApproachInRPhi * closestApproach(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
bool isGoodSecondaryTrack(const reco::TrackRef &primTrack, const reco::TrackRef &secTrack) const
GlobalPoint crossingPoint() const override
bool NuclearVertexBuilder::isGoodSecondaryTrack ( const reco::TrackRef secTrack,
const reco::TrackRef primTrack,
const double &  distOfClosestApp,
const GlobalPoint crossPoint 
) const
private

Definition at line 139 of file NuclearVertexBuilder.cc.

References chi2Cut_, gather_cfg::cout, DPtovPtCut_, minDistFromPrim_, p1, p2, PV3DBase< T, PVType, FrameType >::perp(), and hiDetachedQuadStep_cff::pt2.

142  {
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 }
T perp() const
Definition: PV3DBase.h:72
double p2[4]
Definition: TauolaWrapper.h:90
double p1[4]
Definition: TauolaWrapper.h:89

Member Data Documentation

double NuclearVertexBuilder::chi2Cut_
private

Definition at line 56 of file NuclearVertexBuilder.h.

Referenced by isGoodSecondaryTrack().

double NuclearVertexBuilder::DPtovPtCut_
private

Definition at line 57 of file NuclearVertexBuilder.h.

Referenced by isGoodSecondaryTrack().

double NuclearVertexBuilder::minDistFromPrim_
private

Definition at line 55 of file NuclearVertexBuilder.h.

Referenced by isGoodSecondaryTrack().

double NuclearVertexBuilder::minDistFromVtx_
private

Definition at line 58 of file NuclearVertexBuilder.h.

Referenced by isCompatible().

double NuclearVertexBuilder::shareFrac_
private

Definition at line 59 of file NuclearVertexBuilder.h.

Referenced by cleanTrackCollection().

reco::Vertex NuclearVertexBuilder::the_vertex
private
const MagneticField* NuclearVertexBuilder::theMagField
private

Definition at line 53 of file NuclearVertexBuilder.h.

Referenced by getTrajectory().

const TransientTrackBuilder* NuclearVertexBuilder::theTransientTrackBuilder
private

Definition at line 54 of file NuclearVertexBuilder.h.

Referenced by FillVertexWithAdaptVtxFitter().