CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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  bool outOnly) :
34  theGeometricSearchTracker(tkLayout), theOutLayerOnlyFlag(outOnly), theEpsilon(-0.01) {
35 
36 }
37 
38 
39 //
40 // return compatible layers for a given trajectory state
41 //
42 vector<const DetLayer*>
44  PropagationDirection dir) const {
45 
46  bool inOut = outward(fts);
47  double eta0 = fts.position().eta();
48 
49  vector<const DetLayer*> output;
50 
51  // check eta of DetLayers for compatibility
52 
53  if (inOut) {
54 
55  if ( !theOutLayerOnlyFlag ) {
56  inOutPx(fts,output);
57 
58  if ( eta0 > 1.55) inOutFPx(fts,output);
59  else if ( eta0 < -1.55 ) inOutBPx(fts,output);
60 
61  if ( fabs(eta0) < 1.67 ) inOutTIB(fts,output);
62 
63  if ( eta0 > 1.17 ) inOutFTID(fts,output);
64  else if ( eta0 < -1.17 ) inOutBTID(fts,output);
65  }
66 
67  if ( fabs(eta0) < 1.35 ) inOutTOB(fts,output);
68 
69  if ( eta0 > 0.97 ) inOutFTEC(fts,output);
70  else if ( eta0 < -0.97 ) inOutBTEC(fts,output);
71 
72  } else {
73  LogTrace("Muon|RecoMuon|DirectionTrackerNavigation")<<"No implementation for inward state at this moment. ";
74  }
75 
76  if ( dir == oppositeToMomentum ) std::reverse(output.begin(),output.end());
77 
78  return output;
79 
80 }
81 
82 
83 //
84 //
85 //
86 void DirectTrackerNavigation::inOutPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
87 
88  const vector<BarrelDetLayer*>& barrel = theGeometricSearchTracker->pixelBarrelLayers();
89 
90  for (vector<BarrelDetLayer*>::const_iterator iter_B = barrel.begin(); iter_B != barrel.end(); iter_B++) {
91  if ( checkCompatible(fts,(*iter_B)) ) output.push_back((*iter_B));
92  }
93 
94 }
95 
96 
97 //
98 //
99 //
100 void DirectTrackerNavigation::inOutTIB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
101 
102  const vector<BarrelDetLayer*>& barrel = theGeometricSearchTracker->tibLayers();
103 
104  for (vector<BarrelDetLayer*>::const_iterator iter_B = barrel.begin(); iter_B != barrel.end(); iter_B++) {
105  if ( checkCompatible(fts,(*iter_B)) ) output.push_back((*iter_B));
106  }
107 
108 }
109 
110 
111 //
112 //
113 //
114 void DirectTrackerNavigation::inOutTOB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
115 
116  const vector<BarrelDetLayer*>& barrel = theGeometricSearchTracker->tobLayers();
117 
118  for (vector<BarrelDetLayer*>::const_iterator iter_B = barrel.begin(); iter_B != barrel.end(); iter_B++) {
119  if ( checkCompatible(fts,(*iter_B)) ) output.push_back((*iter_B));
120  }
121 
122 }
123 
124 
125 //
126 //
127 //
128 void DirectTrackerNavigation::inOutFPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
129 
130  const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->posPixelForwardLayers();
131 
132  for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
133  if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
134  }
135 
136 }
137 
138 
139 //
140 //
141 //
142 void DirectTrackerNavigation::inOutFTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
143 
144  const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->posTidLayers();
145 
146  for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
147  if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
148  }
149 
150 }
151 
152 
153 //
154 //
155 //
156 void DirectTrackerNavigation::inOutFTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
157 
158  const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->posTecLayers();
159 
160  for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
161  if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
162  }
163 
164 }
165 
166 
167 //
168 //
169 //
170 void DirectTrackerNavigation::inOutBPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
171 
172  const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->negPixelForwardLayers();
173 
174  for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
175  if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
176  }
177 
178 }
179 
180 
181 //
182 //
183 //
184 void DirectTrackerNavigation::inOutBTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
185 
186  const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->negTidLayers();
187 
188  for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
189  if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
190  }
191 
192 }
193 
194 
195 //
196 //
197 //
198 void DirectTrackerNavigation::inOutBTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
199 
200  const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->negTecLayers();
201 
202  for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
203  if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
204  }
205 
206 }
207 
208 
209 //
210 //
211 //
213 
214  float eta0 = fts.position().eta();
215 
216  const BoundCylinder& bc = dl->specificSurface();
217  float radius = bc.radius();
218  float length = bc.bounds().length()/2.;
219 
220  float eta = calculateEta(radius, length);
221 
222  return ( fabs(eta0) <= (fabs(eta) + theEpsilon) );
223 
224 }
225 
226 
227 //
228 //
229 //
231 
232  float eta0 = fts.position().eta();
233 
234  const BoundDisk& bd = dl->specificSurface();
235 
236  float outRadius = bd.outerRadius();
237  float inRadius = bd.innerRadius();
238  float z = bd.position().z();
239 
240  float etaOut = calculateEta(outRadius, z);
241  float etaIn = calculateEta(inRadius, z);
242 
243  if ( eta0 > 0 ) return ( eta0 > ( etaOut - theEpsilon) && eta0 < (etaIn + theEpsilon) );
244  else return ( eta0 < (etaOut + theEpsilon ) && eta0 > ( etaIn - theEpsilon ) );
245 
246 }
247 
248 
249 //
250 //
251 //
253 
254  return (fts.position().basicVector().dot(fts.momentum().basicVector()) > 0);
255 
256 }
257 
258 
259 //
260 // calculate pseudorapidity from r and z
261 //
262 float DirectTrackerNavigation::calculateEta(float r, float z) const {
263 
264  if ( z > 0 ) return -log((tan(atan(r/z)/2.)));
265  return log(-(tan(atan(r/z)/2.)));
266 
267 }
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
T eta() const
float float float z
void inOutTOB(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
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
#define LogTrace(id)
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
bool outward(const FreeTrajectoryState &) const
void inOutFTEC(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
virtual const BoundDisk & specificSurface() const GCC11_FINAL
T eta() const
Definition: PV3DBase.h:76
std::vector< const DetLayer * > compatibleLayers(const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
find compatible layers for a given trajectory state
dbl *** dir
Definition: mlp_gen.cc:35
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
virtual const BoundCylinder & specificSurface() const GCC11_FINAL
Extension of the interface.
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.