test
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  for (const auto i : theGeometricSearchTracker->pixelBarrelLayers()) {
89  if ( checkCompatible(fts,i) ) output.push_back(i);
90  }
91 
92 }
93 
94 
95 //
96 //
97 //
98 void DirectTrackerNavigation::inOutTIB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
99 
100  for(const auto i : theGeometricSearchTracker->tibLayers() ) {
101  if ( checkCompatible(fts,i) ) output.push_back(i);
102  }
103 
104 }
105 
106 
107 //
108 //
109 //
110 void DirectTrackerNavigation::inOutTOB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
111 
112  for( const auto i : theGeometricSearchTracker->tobLayers()) {
113  if ( checkCompatible(fts,i) ) output.push_back(i);
114  }
115 
116 }
117 
118 
119 //
120 //
121 //
122 void DirectTrackerNavigation::inOutFPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
123 
124  for(const auto i : theGeometricSearchTracker->posPixelForwardLayers()) {
125  if ( checkCompatible(fts,i) ) output.push_back(i);
126  }
127 
128 }
129 
130 
131 //
132 //
133 //
134 void DirectTrackerNavigation::inOutFTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
135 
136  for(const auto i: theGeometricSearchTracker->posTidLayers()) {
137  if ( checkCompatible(fts,i) ) output.push_back(i);
138  }
139 
140 }
141 
142 
143 //
144 //
145 //
146 void DirectTrackerNavigation::inOutFTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
147 
148  for(const auto i : theGeometricSearchTracker->posTecLayers()) {
149  if ( checkCompatible(fts,i) ) output.push_back(i);
150  }
151 
152 }
153 
154 
155 //
156 //
157 //
158 void DirectTrackerNavigation::inOutBPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
159 
160  for(const auto i: theGeometricSearchTracker->negPixelForwardLayers()) {
161  if ( checkCompatible(fts,i) ) output.push_back(i);
162  }
163 
164 }
165 
166 
167 //
168 //
169 //
170 void DirectTrackerNavigation::inOutBTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
171 
172  for(const auto i: theGeometricSearchTracker->negTidLayers()) {
173  if ( checkCompatible(fts,i) ) output.push_back(i);
174  }
175 
176 }
177 
178 
179 //
180 //
181 //
182 void DirectTrackerNavigation::inOutBTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
183 
184  for(const auto i: theGeometricSearchTracker->negTecLayers()) {
185  if ( checkCompatible(fts,i) ) output.push_back(i);
186  }
187 
188 }
189 
190 
191 //
192 //
193 //
195 
196  float eta0 = fts.position().eta();
197 
198  const BoundCylinder& bc = dl->specificSurface();
199  float radius = bc.radius();
200  float length = bc.bounds().length()/2.;
201 
202  float eta = calculateEta(radius, length);
203 
204  return ( fabs(eta0) <= (fabs(eta) + theEpsilon) );
205 
206 }
207 
208 
209 //
210 //
211 //
213 
214  float eta0 = fts.position().eta();
215 
216  const BoundDisk& bd = dl->specificSurface();
217 
218  float outRadius = bd.outerRadius();
219  float inRadius = bd.innerRadius();
220  float z = bd.position().z();
221 
222  float etaOut = calculateEta(outRadius, z);
223  float etaIn = calculateEta(inRadius, z);
224 
225  if ( eta0 > 0 ) return ( eta0 > ( etaOut - theEpsilon) && eta0 < (etaIn + theEpsilon) );
226  else return ( eta0 < (etaOut + theEpsilon ) && eta0 > ( etaIn - theEpsilon ) );
227 
228 }
229 
230 
231 //
232 //
233 //
235 
236  return (fts.position().basicVector().dot(fts.momentum().basicVector()) > 0);
237 
238 }
239 
240 
241 //
242 // calculate pseudorapidity from r and z
243 //
244 float DirectTrackerNavigation::calculateEta(float r, float z) const {
245 
246  if ( z > 0 ) return -log((tan(atan(r/z)/2.)));
247  return log(-(tan(atan(r/z)/2.)));
248 
249 }
int i
Definition: DBlmapReader.cc:9
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
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
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:76
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
dbl *** dir
Definition: mlp_gen.cc:35
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
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.