CMS 3D CMS Logo

FWPFTrackUtils.cc
Go to the documentation of this file.
2 
5 
6 //______________________________________________________________________________
8  if (!instanceFlag) // Instance doesn't exist yet
9  {
11  instanceFlag = true;
12  }
13 
14  return pInstance; // Pointer to sole instance
15 }
16 
17 //______________________________________________________________________________
19  m_magField = new FWMagField();
20 
21  // Common propagator, helix stepper
22  m_trackPropagator = new TEveTrackPropagator();
23  m_trackPropagator->SetMagFieldObj(m_magField, false);
26  m_trackPropagator->SetDelta(0.01);
27  m_trackPropagator->SetProjTrackBreaking(TEveTrackPropagator::kPTB_UseLastPointPos);
28  m_trackPropagator->SetRnrPTBMarkers(kTRUE);
29  m_trackPropagator->IncDenyDestroy();
30 
31  // Tracker propagator
32  m_trackerTrackPropagator = new TEveTrackPropagator();
33  m_trackerTrackPropagator->SetStepper(TEveTrackPropagator::kRungeKutta);
34  m_trackerTrackPropagator->SetMagFieldObj(m_magField, false);
35  m_trackerTrackPropagator->SetDelta(0.01);
38  m_trackerTrackPropagator->SetProjTrackBreaking(TEveTrackPropagator::kPTB_UseLastPointPos);
39  m_trackerTrackPropagator->SetRnrPTBMarkers(kTRUE);
40  m_trackerTrackPropagator->IncDenyDestroy();
41 }
42 
43 //______________________________________________________________________________
45 
46 //______________________________________________________________________________
47 TEveTrack *FWPFTrackUtils::getTrack(const reco::Track &iData) {
48  TEveTrackPropagator *propagator =
50 
51  TEveRecTrack t;
52  t.fBeta = 1;
53  t.fP = TEveVector(iData.px(), iData.py(), iData.pz());
54  t.fV = TEveVector(iData.vertex().x(), iData.vertex().y(), iData.vertex().z());
55  t.fSign = iData.charge();
56  TEveTrack *trk = new TEveTrack(&t, propagator);
57  trk->MakeTrack();
58 
59  return trk;
60 }
61 
62 //______________________________________________________________________________
63 TEveStraightLineSet *FWPFTrackUtils::setupLegoTrack(const reco::Track &iData) {
64  using namespace FWPFGeom;
65 
66  // Declarations
67  int wraps[3] = {-1, -1, -1};
68  bool ECAL = false;
69  TEveTrack *trk = getTrack(iData);
70  std::vector<TEveVector> trackPoints(trk->GetN() - 1);
71  const Float_t *points = trk->GetP();
72  TEveStraightLineSet *legoTrack = new TEveStraightLineSet();
73 
75  if (fabs(iData.eta()) < 2.0 && iData.pt() > 0.5 && iData.pt() < 30) {
76  double estimate = fw::estimate_field(iData, true);
77  if (estimate >= 0)
78  m_singleton->getField()->guessField(estimate);
79  }
80  }
81 
82  // Convert to Eta/Phi and store in vector
83  for (Int_t i = 1; i < trk->GetN(); ++i) {
84  int j = i * 3;
85  TEveVector temp = TEveVector(points[j], points[j + 1], points[j + 2]);
86  TEveVector vec = TEveVector(temp.Eta(), temp.Phi(), 0.001);
87 
88  trackPoints[i - 1] = vec;
89  }
90 
91  // Add first point to ps if necessary
92  for (Int_t i = 1; i < trk->GetN(); ++i) {
93  int j = i * 3;
94  TEveVector v1 = TEveVector(points[j], points[j + 1], points[j + 2]);
95 
96  if (!ECAL) {
97  if (FWPFMaths::checkIntersect(v1, caloR1())) {
98  TEveVector v2 = TEveVector(points[j - 3], points[j - 2], points[j - 1]);
99  TEveVector xyPoint = FWPFMaths::lineCircleIntersect(v1, v2, caloR1());
100  TEveVector zPoint;
101  if (v1.fZ < 0)
102  zPoint = TEveVector(xyPoint.fX, xyPoint.fY, v1.fZ - 50.f);
103  else
104  zPoint = TEveVector(xyPoint.fX, xyPoint.fY, v1.fZ + 50.f);
105 
106  TEveVector vec = FWPFMaths::lineLineIntersect(v1, v2, xyPoint, zPoint);
107  legoTrack->AddMarker(vec.Eta(), vec.Phi(), 0.001, 0);
108 
109  wraps[0] = i; // There is now a chance that the track will also reach the HCAL radius
110  ECAL = true;
111  } else if (fabs(v1.fZ) >= caloZ1()) {
112  TEveVector p1, p2;
113  TEveVector vec, v2 = TEveVector(points[j - 3], points[j - 2], points[j - 1]);
114  float z, y = FWPFMaths::linearInterpolation(v2, v1, caloZ1());
115 
116  if (v2.fZ > 0)
117  z = caloZ1();
118  else
119  z = caloZ1() * -1;
120 
121  p1 = TEveVector(v2.fX + 50, y, z);
122  p2 = TEveVector(v2.fX - 50, y, z);
123  vec = FWPFMaths::lineLineIntersect(v1, v2, p1, p2);
124 
125  legoTrack->AddMarker(vec.Eta(), vec.Phi(), 0.001, 0);
126  wraps[0] = i;
127  ECAL = true;
128  }
129  } else if (FWPFMaths::checkIntersect(v1, caloR2())) {
130  TEveVector v2 = TEveVector(points[j - 3], points[j - 2], points[j - 1]);
131  TEveVector xyPoint = FWPFMaths::lineCircleIntersect(v1, v2, caloR2());
132  TEveVector zPoint;
133  if (v1.fZ < 0)
134  zPoint = TEveVector(xyPoint.fX, xyPoint.fY, v1.fZ - 50.f);
135  else
136  zPoint = TEveVector(xyPoint.fX, xyPoint.fY, v1.fZ + 50.f);
137 
138  TEveVector vec = FWPFMaths::lineLineIntersect(v1, v2, xyPoint, zPoint);
139  legoTrack->AddMarker(vec.Eta(), vec.Phi(), 0.001, 1);
140 
141  wraps[1] = i; // There is now a chance that the track will also reach the HCAL radius
142  break;
143  }
144  }
145 
146  if (wraps[0] != -1) //if( ECAL )
147  {
148  int i = (trk->GetN() - 1) * 3;
149  int j = trk->GetN() - 2; // This is equal to the last index in trackPoints vector
150  TEveVector vec = TEveVector(points[i], points[i + 1], points[i + 2]);
151 
152  if (FWPFMaths::checkIntersect(vec, caloR3() - 1)) {
153  legoTrack->AddMarker(vec.Eta(), vec.Phi(), 0.001, 2);
154  wraps[2] = j;
155  } else if (fabs(vec.fZ) >= caloZ2()) {
156  legoTrack->AddMarker(vec.Eta(), vec.Phi(), 0.001, 2);
157  wraps[2] = j;
158  }
159  }
160 
161  /* Handle phi wrapping */
162  for (int i = 0; i < static_cast<int>(trackPoints.size() - 1); ++i) {
163  if ((trackPoints[i + 1].fY - trackPoints[i].fY) > 1) {
164  trackPoints[i + 1].fY -= TMath::TwoPi();
165  if (i == wraps[0]) {
166  TEveChunkManager::iterator mi(legoTrack->GetMarkerPlex());
167  mi.next(); // First point
168  TEveStraightLineSet::Marker_t &m = *(TEveStraightLineSet::Marker_t *)mi();
169  m.fV[0] = trackPoints[i + 1].fX;
170  m.fV[1] = trackPoints[i + 1].fY;
171  m.fV[2] = 0.001;
172  } else if (i == wraps[1]) {
173  TEveChunkManager::iterator mi(legoTrack->GetMarkerPlex());
174  mi.next();
175  mi.next(); // Second point
176  TEveStraightLineSet::Marker_t &m = *(TEveStraightLineSet::Marker_t *)mi();
177  m.fV[0] = trackPoints[i + 1].fX;
178  m.fV[1] = trackPoints[i + 1].fY;
179  m.fV[2] = 0.001;
180  }
181  }
182 
183  if ((trackPoints[i].fY - trackPoints[i + 1].fY) > 1) {
184  trackPoints[i + 1].fY += TMath::TwoPi();
185  if (i == wraps[0]) {
186  TEveChunkManager::iterator mi(legoTrack->GetMarkerPlex());
187  mi.next(); // First point
188  TEveStraightLineSet::Marker_t &m = *(TEveStraightLineSet::Marker_t *)mi();
189  m.fV[0] = trackPoints[i + 1].fX;
190  m.fV[1] = trackPoints[i + 1].fY;
191  m.fV[2] = 0.001;
192  } else if (i == wraps[1]) {
193  TEveChunkManager::iterator mi(legoTrack->GetMarkerPlex());
194  mi.next();
195  mi.next(); // Second point
196  TEveStraightLineSet::Marker_t &m = *(TEveStraightLineSet::Marker_t *)mi();
197  m.fV[0] = trackPoints[i + 1].fX;
198  m.fV[1] = trackPoints[i + 1].fY;
199  m.fV[2] = 0.001;
200  }
201  }
202  }
203 
204  int end = static_cast<int>(trackPoints.size() - 1);
205  if (wraps[2] == end) {
206  TEveChunkManager::iterator mi(legoTrack->GetMarkerPlex());
207  mi.next();
208  mi.next();
209  mi.next(); // Third point
210  TEveStraightLineSet::Marker_t &m = *(TEveStraightLineSet::Marker_t *)mi();
211  m.fV[0] = trackPoints[end].fX;
212  m.fV[1] = trackPoints[end].fY;
213  m.fV[2] = 0.001;
214  }
215 
216  // Set points on TEveLineSet object ready for displaying
217  for (unsigned int i = 0; i < trackPoints.size() - 1; ++i)
218  legoTrack->AddLine(trackPoints[i], trackPoints[i + 1]);
219 
220  legoTrack->SetDepthTest(false);
221  legoTrack->SetMarkerStyle(4);
222  legoTrack->SetMarkerSize(1);
223  legoTrack->SetRnrMarkers(true);
224 
225  delete trk; // Release memory that is no longer required
226 
227  return legoTrack;
228 }
229 
230 //______________________________________________________________________________
231 TEveTrack *FWPFTrackUtils::setupTrack(const reco::Track &iData) {
233  if (fabs(iData.eta()) < 2.0 && iData.pt() > 0.5 && iData.pt() < 30) {
234  double estimate = fw::estimate_field(iData, true);
235  if (estimate >= 0)
236  m_singleton->getField()->guessField(estimate);
237  }
238  }
239 
240  TEveTrack *trk = getTrack(iData);
241 
242  return trk;
243 }
244 
245 //______________________________________________________________________________
246 TEvePointSet *FWPFTrackUtils::getCollisionMarkers(const TEveTrack *trk) {
247  using namespace FWPFGeom;
248 
249  bool ECAL = false;
250  const Float_t *points = trk->GetP();
251  TEvePointSet *ps = new TEvePointSet();
252 
253  for (Int_t i = 1; i < trk->GetN(); ++i) {
254  int j = i * 3;
255  TEveVector v1 = TEveVector(points[j], points[j + 1], points[j + 2]);
256 
257  if (!ECAL) {
258  if (FWPFMaths::checkIntersect(v1, caloR1())) {
259  TEveVector v2 = TEveVector(points[j - 3], points[j - 2], points[j - 1]);
260  TEveVector xyPoint = FWPFMaths::lineCircleIntersect(v1, v2, caloR1());
261  TEveVector zPoint;
262  if (v1.fZ < 0)
263  zPoint = TEveVector(xyPoint.fX, xyPoint.fY, v1.fZ - 50.f);
264  else
265  zPoint = TEveVector(xyPoint.fX, xyPoint.fY, v1.fZ + 50.f);
266 
267  TEveVector vec = FWPFMaths::lineLineIntersect(v1, v2, xyPoint, zPoint);
268  ps->SetNextPoint(vec.fX, vec.fY, vec.fZ);
269 
270  ECAL = true;
271  } else if (fabs(v1.fZ) >= caloZ1()) {
272  TEveVector p1, p2;
273  TEveVector vec, v2 = TEveVector(points[j - 3], points[j - 2], points[j - 1]);
274  float z, y = FWPFMaths::linearInterpolation(v2, v1, caloZ1());
275 
276  if (v2.fZ > 0)
277  z = caloZ1();
278  else
279  z = caloZ1() * -1;
280 
281  p1 = TEveVector(v2.fX + 50, y, z);
282  p2 = TEveVector(v2.fX - 50, y, z);
283  vec = FWPFMaths::lineLineIntersect(v1, v2, p1, p2);
284 
285  ps->SetNextPoint(vec.fX, vec.fY, vec.fZ);
286  ECAL = true;
287  }
288  } else if (FWPFMaths::checkIntersect(v1, caloR2())) {
289  TEveVector v2 = TEveVector(points[j - 3], points[j - 2], points[j - 1]);
290  TEveVector xyPoint = FWPFMaths::lineCircleIntersect(v1, v2, caloR2());
291  TEveVector zPoint;
292  if (v1.fZ < 0)
293  zPoint = TEveVector(xyPoint.fX, xyPoint.fY, v1.fZ - 50.f);
294  else
295  zPoint = TEveVector(xyPoint.fX, xyPoint.fY, v1.fZ + 50.f);
296 
297  TEveVector vec = FWPFMaths::lineLineIntersect(v1, v2, xyPoint, zPoint);
298  ps->SetNextPoint(vec.fX, vec.fY, vec.fZ);
299  break; // ECAL and HCAL collisions found so stop looping
300  }
301  }
302 
303  // HCAL collisions (outer radius and endcap)
304  int i = (trk->GetN() - 1) * 3;
305  TEveVector vec = TEveVector(points[i], points[i + 1], points[i + 2]);
306 
307  if (FWPFMaths::checkIntersect(vec, caloR3() - 1))
308  ps->SetNextPoint(vec.fX, vec.fY, vec.fZ);
309  else if (fabs(vec.fZ) >= caloZ2())
310  ps->SetNextPoint(vec.fX, vec.fY, vec.fZ);
311 
312  return ps;
313 }
FWPFTrackUtils.h
DDAxes::y
FWPFTrackUtils::m_singleton
FWPFTrackSingleton * m_singleton
Definition: FWPFTrackUtils.h:85
FWMagField::guessField
void guessField(float estimate) const
Definition: FWMagField.cc:129
mps_fire.i
i
Definition: mps_fire.py:428
FWPFTrackUtils::FWPFTrackUtils
FWPFTrackUtils()
Definition: FWPFTrackUtils.cc:44
TwoPi
const double TwoPi
Definition: CosmicMuonParameters.h:19
FWPFTrackSingleton::getField
FWMagField * getField()
Definition: FWPFTrackUtils.h:41
FWMagField::kNone
Definition: FWMagField.h:18
FWPFTrackSingleton::m_trackPropagator
TEveTrackPropagator * m_trackPropagator
Definition: FWPFTrackUtils.h:60
HLT_FULL_cff.points
points
Definition: HLT_FULL_cff.py:21449
edm::Ref::isAvailable
bool isAvailable() const
Definition: Ref.h:537
FWPFTrackUtils::setupLegoTrack
TEveStraightLineSet * setupLegoTrack(const reco::Track &)
Definition: FWPFTrackUtils.cc:63
FWPFTrackUtils::getTrack
TEveTrack * getTrack(const reco::Track &)
Definition: FWPFTrackUtils.cc:47
FWPFTrackSingleton::Instance
static FWPFTrackSingleton * Instance()
Definition: FWPFTrackUtils.cc:7
reco::TrackBase::px
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:640
FWPFTrackSingleton::initPropagator
void initPropagator()
Definition: FWPFTrackUtils.cc:18
FWPFTrackSingleton::getTrackerTrackPropagator
TEveTrackPropagator * getTrackerTrackPropagator()
Definition: FWPFTrackUtils.h:39
groupFilesInBlocks.temp
list temp
Definition: groupFilesInBlocks.py:142
reco::Track::extra
const TrackExtraRef & extra() const
reference to "extra" object
Definition: Track.h:139
FWPFMaths::lineLineIntersect
TEveVector lineLineIntersect(const TEveVector &v1, const TEveVector &v2, const TEveVector &v3, const TEveVector &v4)
Definition: FWPFMaths.cc:73
FWPFTrackUtils::getCollisionMarkers
TEvePointSet * getCollisionMarkers(const TEveTrack *)
Definition: FWPFTrackUtils.cc:246
reco::TrackBase::pt
double pt() const
track transverse momentum
Definition: TrackBase.h:637
ECAL
Definition: HCALResponse.h:21
TrackCandidateProducer_cfi.propagator
propagator
Definition: TrackCandidateProducer_cfi.py:17
reco::TrackBase::py
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:643
FWPFGeom::caloR3
float caloR3()
Definition: FWPFGeom.h:27
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
mps_fire.end
end
Definition: mps_fire.py:242
FWMagField
Definition: FWMagField.h:14
DDAxes::z
reco::Track
Definition: Track.h:27
p2
double p2[4]
Definition: TauolaWrapper.h:90
FWMagField::getSource
ESource getSource() const
Definition: FWMagField.h:31
reco::TrackBase::charge
int charge() const
track electric charge
Definition: TrackBase.h:596
FWPFTrackSingleton::m_magField
FWMagField * m_magField
Definition: FWPFTrackUtils.h:61
FWPFGeom::caloR1
float caloR1()
Definition: FWPFGeom.h:20
FWPFTrackUtils::setupTrack
TEveTrack * setupTrack(const reco::Track &)
Definition: FWPFTrackUtils.cc:231
reco::TrackBase::eta
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
FWPFGeom::caloZ1
float caloZ1()
Definition: FWPFGeom.h:21
FWPFMaths::checkIntersect
bool checkIntersect(const TEveVector &vec, float r)
Definition: FWPFMaths.cc:103
p1
double p1[4]
Definition: TauolaWrapper.h:89
FWPFTrackSingleton::getTrackPropagator
TEveTrackPropagator * getTrackPropagator()
Definition: FWPFTrackUtils.h:40
reco::TrackBase::vertex
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead
Definition: TrackBase.h:676
FWPFMaths::lineCircleIntersect
TEveVector lineCircleIntersect(const TEveVector &v1, const TEveVector &v2, float r)
Definition: FWPFMaths.cc:38
FWPFGeom::caloZ2
float caloZ2()
Definition: FWPFGeom.h:28
FWPFTrackSingleton::instanceFlag
static bool instanceFlag
Definition: FWPFTrackUtils.h:57
FWPFTrackSingleton::pInstance
static FWPFTrackSingleton * pInstance
Definition: FWPFTrackUtils.h:56
FWPFMaths::linearInterpolation
float linearInterpolation(const TEveVector &p1, const TEveVector &p2, float r)
Definition: FWPFMaths.cc:91
reco::TrackBase::pz
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:646
FWPFGeom
Definition: FWPFGeom.h:18
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
FWPFTrackSingleton::FWPFTrackSingleton
FWPFTrackSingleton()
Definition: FWPFTrackUtils.h:48
fw::estimate_field
double estimate_field(const reco::Track &track, bool highQuality=false)
Definition: estimate_field.cc:11
FWPFGeom::caloR2
float caloR2()
Definition: FWPFGeom.h:26
FWPFTrackSingleton::m_trackerTrackPropagator
TEveTrackPropagator * m_trackerTrackPropagator
Definition: FWPFTrackUtils.h:59
FWPFTrackSingleton
Definition: FWPFTrackUtils.h:34