CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DirectTrackerNavigation.cc
Go to the documentation of this file.
1 
10 
11 //---------------
12 // C++ Headers --
13 //---------------
14 
15 #include <algorithm>
16 
17 //-------------------------------
18 // Collaborating Class Headers --
19 //-------------------------------
20 
26 
27 using namespace std;
28 
29 //
30 // constructor
31 //
33  : theGeometricSearchTracker(tkLayout), theOutLayerOnlyFlag(outOnly), theEpsilon(-0.01) {}
34 
35 //
36 // return compatible layers for a given trajectory state
37 //
39  PropagationDirection dir) const {
40  bool inOut = outward(fts);
41  double eta0 = fts.position().eta();
42 
43  vector<const DetLayer*> output;
44 
45  // check eta of DetLayers for compatibility
46 
47  if (inOut) {
48  if (!theOutLayerOnlyFlag) {
49  inOutPx(fts, output);
50 
51  if (eta0 > 1.55)
52  inOutFPx(fts, output);
53  else if (eta0 < -1.55)
54  inOutBPx(fts, output);
55 
56  if (fabs(eta0) < 1.67)
57  inOutTIB(fts, output);
58 
59  if (eta0 > 1.17)
60  inOutFTID(fts, output);
61  else if (eta0 < -1.17)
62  inOutBTID(fts, output);
63  }
64 
65  if (fabs(eta0) < 1.35)
66  inOutTOB(fts, output);
67 
68  if (eta0 > 0.97)
69  inOutFTEC(fts, output);
70  else if (eta0 < -0.97)
71  inOutBTEC(fts, output);
72 
73  } else {
74  LogTrace("Muon|RecoMuon|DirectionTrackerNavigation") << "No implementation for inward state at this moment. ";
75  }
76 
77  if (dir == oppositeToMomentum)
78  std::reverse(output.begin(), output.end());
79 
80  return output;
81 }
82 
83 //
84 //
85 //
86 void DirectTrackerNavigation::inOutPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
87  for (const auto i : theGeometricSearchTracker->pixelBarrelLayers()) {
88  if (checkCompatible(fts, i))
89  output.push_back(i);
90  }
91 }
92 
93 //
94 //
95 //
96 void DirectTrackerNavigation::inOutTIB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
97  for (const auto i : theGeometricSearchTracker->tibLayers()) {
98  if (checkCompatible(fts, i))
99  output.push_back(i);
100  }
101 }
102 
103 //
104 //
105 //
106 void DirectTrackerNavigation::inOutTOB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
107  for (const auto i : theGeometricSearchTracker->tobLayers()) {
108  if (checkCompatible(fts, i))
109  output.push_back(i);
110  }
111 }
112 
113 //
114 //
115 //
116 void DirectTrackerNavigation::inOutFPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
117  for (const auto i : theGeometricSearchTracker->posPixelForwardLayers()) {
118  if (checkCompatible(fts, i))
119  output.push_back(i);
120  }
121 }
122 
123 //
124 //
125 //
126 void DirectTrackerNavigation::inOutFTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
127  for (const auto i : theGeometricSearchTracker->posTidLayers()) {
128  if (checkCompatible(fts, i))
129  output.push_back(i);
130  }
131 }
132 
133 //
134 //
135 //
136 void DirectTrackerNavigation::inOutFTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
137  for (const auto i : theGeometricSearchTracker->posTecLayers()) {
138  if (checkCompatible(fts, i))
139  output.push_back(i);
140  }
141 }
142 
143 //
144 //
145 //
146 void DirectTrackerNavigation::inOutBPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
147  for (const auto i : theGeometricSearchTracker->negPixelForwardLayers()) {
148  if (checkCompatible(fts, i))
149  output.push_back(i);
150  }
151 }
152 
153 //
154 //
155 //
156 void DirectTrackerNavigation::inOutBTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
157  for (const auto i : theGeometricSearchTracker->negTidLayers()) {
158  if (checkCompatible(fts, i))
159  output.push_back(i);
160  }
161 }
162 
163 //
164 //
165 //
166 void DirectTrackerNavigation::inOutBTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
167  for (const auto i : theGeometricSearchTracker->negTecLayers()) {
168  if (checkCompatible(fts, i))
169  output.push_back(i);
170  }
171 }
172 
173 //
174 //
175 //
177  float eta0 = fts.position().eta();
178 
179  const BoundCylinder& bc = dl->specificSurface();
180  float radius = bc.radius();
181  float length = bc.bounds().length() / 2.;
182 
183  float eta = calculateEta(radius, length);
184 
185  return (fabs(eta0) <= (fabs(eta) + theEpsilon));
186 }
187 
188 //
189 //
190 //
192  float eta0 = fts.position().eta();
193 
194  const BoundDisk& bd = dl->specificSurface();
195 
196  float outRadius = bd.outerRadius();
197  float inRadius = bd.innerRadius();
198  float z = bd.position().z();
199 
200  float etaOut = calculateEta(outRadius, z);
201  float etaIn = calculateEta(inRadius, z);
202 
203  if (eta0 > 0)
204  return (eta0 > (etaOut - theEpsilon) && eta0 < (etaIn + theEpsilon));
205  else
206  return (eta0 < (etaOut + theEpsilon) && eta0 > (etaIn - theEpsilon));
207 }
208 
209 //
210 //
211 //
213  return (fts.position().basicVector().dot(fts.momentum().basicVector()) > 0);
214 }
215 
216 //
217 // calculate pseudorapidity from r and z
218 //
219 float DirectTrackerNavigation::calculateEta(float r, float z) const {
220  if (z > 0)
221  return -log((tan(atan(r / z) / 2.)));
222  return log(-(tan(atan(r / z) / 2.)));
223 }
static std::vector< std::string > checklist log
void inOutPx(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
void inOutTIB(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
void inOutFTID(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
void inOutBTID(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
PropagationDirection
float float float z
void inOutTOB(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
#define LogTrace(id)
edm::ESHandle< GeometricSearchTracker > theGeometricSearchTracker
DirectTrackerNavigation(const edm::ESHandle< GeometricSearchTracker > &, bool outOnly=true)
constructor
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
bool checkCompatible(const FreeTrajectoryState &, const BarrelDetLayer *) const
GlobalVector momentum() const
void inOutFPx(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
GlobalPoint position() const
float calculateEta(float r, float z) const
void inOutBTEC(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.
bool outward(const FreeTrajectoryState &) const
void inOutFTEC(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
T eta() const
Definition: PV3DBase.h:73
std::vector< const DetLayer * > compatibleLayers(const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
find compatible layers for a given trajectory state
virtual const BoundDisk & specificSurface() const final
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:53
void inOutBPx(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
T dot(const Basic3DVector &rh) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.