CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DirectMuonNavigation.cc
Go to the documentation of this file.
2 
14 
15 #include <algorithm>
16 
17 using namespace std;
18 
20  : theMuonDetLayerGeometry(muonLayout), epsilon_(100.), theEndcapFlag(true), theBarrelFlag(true) {}
21 
23  const edm::ParameterSet& par)
24  : theMuonDetLayerGeometry(muonLayout),
25  epsilon_(100.),
26  theEndcapFlag(par.getParameter<bool>("Endcap")),
27  theBarrelFlag(par.getParameter<bool>("Barrel")) {}
28 
29 /* return compatible layers for given trajectory state */
31  PropagationDirection dir) const {
32  float z0 = fts.position().z();
33  float zm = fts.momentum().z();
34 
35  bool inOut = outward(fts);
36 
37  vector<const DetLayer*> output;
38 
39  // check direction and position of FTS to get a correct order of DetLayers
40 
41  if (inOut) {
42  if ((zm * z0) >= 0) {
43  if (theBarrelFlag)
44  inOutBarrel(fts, output);
45  if (theEndcapFlag) {
46  if (z0 >= 0)
47  inOutForward(fts, output);
48  else
49  inOutBackward(fts, output);
50  }
51  } else {
52  if (theEndcapFlag) {
53  if (z0 >= 0)
54  outInForward(fts, output);
55  else
56  outInBackward(fts, output);
57  }
58  if (theBarrelFlag)
59  inOutBarrel(fts, output);
60  }
61  } else {
62  if ((zm * z0) >= 0) {
63  if (theBarrelFlag)
64  outInBarrel(fts, output);
65  if (theEndcapFlag) {
66  if (z0 >= 0)
67  inOutForward(fts, output);
68  else
69  inOutBackward(fts, output);
70  }
71  } else {
72  if (theEndcapFlag) {
73  if (z0 >= 0)
74  outInForward(fts, output);
75  else
76  outInBackward(fts, output);
77  }
78  if (theBarrelFlag)
79  outInBarrel(fts, output);
80  }
81  }
82 
83  if (dir == oppositeToMomentum)
84  std::reverse(output.begin(), output.end());
85 
86  return output;
87 }
88 
89 /*
90 return compatible endcap layers on BOTH ends;
91 used for beam-halo muons
92 */
94  PropagationDirection dir) const {
95  float zm = fts.momentum().z();
96 
97  vector<const DetLayer*> output;
98 
99  // collect all endcap layers on 2 sides
100  outInBackward(fts, output);
101  inOutForward(fts, output);
102 
103  // check direction FTS to get a correct order of DetLayers
104  if ((zm > 0 && dir == oppositeToMomentum) || (zm < 0 && dir == alongMomentum))
105  std::reverse(output.begin(), output.end());
106 
107  return output;
108 }
109 
110 void DirectMuonNavigation::inOutBarrel(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
111  bool cont = false;
112  const vector<const DetLayer*>& barrel = theMuonDetLayerGeometry->allBarrelLayers();
113 
114  for (vector<const DetLayer*>::const_iterator iter_B = barrel.begin(); iter_B != barrel.end(); iter_B++) {
115  if (cont)
116  output.push_back((*iter_B));
117  else if (checkCompatible(fts, dynamic_cast<const BarrelDetLayer*>(*iter_B))) {
118  output.push_back((*iter_B));
119  cont = true;
120  }
121  }
122 }
123 
124 void DirectMuonNavigation::outInBarrel(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
125  // default barrel layers are in out, reverse order
126  const vector<const DetLayer*>& barrel = theMuonDetLayerGeometry->allBarrelLayers();
127 
128  bool cont = false;
129  vector<const DetLayer*>::const_iterator rbegin = barrel.end();
130  rbegin--;
131  vector<const DetLayer*>::const_iterator rend = barrel.begin();
132  rend--;
133 
134  for (vector<const DetLayer*>::const_iterator iter_B = rbegin; iter_B != rend; iter_B--) {
135  if (cont)
136  output.push_back((*iter_B));
137  else if (checkCompatible(fts, dynamic_cast<const BarrelDetLayer*>(*iter_B))) {
138  output.push_back((*iter_B));
139  cont = true;
140  }
141  }
142 }
143 
144 void DirectMuonNavigation::inOutForward(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
145  const vector<const DetLayer*>& forward = theMuonDetLayerGeometry->allForwardLayers();
146  bool cont = false;
147  for (vector<const DetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
148  if (cont)
149  output.push_back((*iter_E));
150  else if (checkCompatible(fts, dynamic_cast<const ForwardDetLayer*>(*iter_E))) {
151  output.push_back((*iter_E));
152  cont = true;
153  }
154  }
155 }
156 
157 void DirectMuonNavigation::outInForward(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
158  // default forward layers are in out, reverse order
159 
160  bool cont = false;
161  const vector<const DetLayer*>& forward = theMuonDetLayerGeometry->allForwardLayers();
162  vector<const DetLayer*>::const_iterator rbegin = forward.end();
163  rbegin--;
164  vector<const DetLayer*>::const_iterator rend = forward.begin();
165  rend--;
166  for (vector<const DetLayer*>::const_iterator iter_E = rbegin; iter_E != rend; iter_E--) {
167  if (cont)
168  output.push_back((*iter_E));
169  else if (checkCompatible(fts, dynamic_cast<const ForwardDetLayer*>(*iter_E))) {
170  output.push_back((*iter_E));
171  cont = true;
172  }
173  }
174 }
175 
176 void DirectMuonNavigation::inOutBackward(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
177  bool cont = false;
178  const vector<const DetLayer*>& backward = theMuonDetLayerGeometry->allBackwardLayers();
179 
180  for (vector<const DetLayer*>::const_iterator iter_E = backward.begin(); iter_E != backward.end(); iter_E++) {
181  if (cont)
182  output.push_back((*iter_E));
183  else if (checkCompatible(fts, dynamic_cast<const ForwardDetLayer*>(*iter_E))) {
184  output.push_back((*iter_E));
185  cont = true;
186  }
187  }
188 }
189 
190 void DirectMuonNavigation::outInBackward(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
191  bool cont = false;
192  const vector<const DetLayer*>& backward = theMuonDetLayerGeometry->allBackwardLayers();
193 
194  vector<const DetLayer*>::const_iterator rbegin = backward.end();
195  rbegin--;
196  vector<const DetLayer*>::const_iterator rend = backward.begin();
197  rend--;
198  for (vector<const DetLayer*>::const_iterator iter_E = rbegin; iter_E != rend; iter_E--) {
199  if (cont)
200  output.push_back((*iter_E));
201  else if (checkCompatible(fts, dynamic_cast<const ForwardDetLayer*>(*iter_E))) {
202  output.push_back((*iter_E));
203  cont = true;
204  }
205  }
206 }
207 
209  float z0 = fts.position().z();
210  float r0 = fts.position().perp();
211  float zm = fts.momentum().z();
212  float rm = fts.momentum().perp();
213  float slope = zm / rm;
214  if (!outward(fts))
215  slope = -slope;
216  const BoundCylinder& bc = dl->specificSurface();
217  float radius = bc.radius();
218  float length = bc.bounds().length() / 2.;
219 
220  float z1 = slope * (radius - r0) + z0;
221  return (fabs(z1) <= fabs(length) + epsilon_);
222 }
223 
225  float z0 = fts.position().z();
226  float r0 = fts.position().perp();
227  float zm = fts.momentum().z();
228  float rm = fts.momentum().perp();
229  float slope = rm / zm;
230 
231  if (!outward(fts))
232  slope = -slope;
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 r1 = slope * (z - z0) + r0;
241  return (r1 >= inRadius - epsilon_ && r1 <= outRadius + epsilon_);
242 }
243 
245  // return (fts.position().basicVector().dot(fts.momentum().basicVector())>0);
246 
247  float x0 = fts.position().x();
248  float y0 = fts.position().y();
249 
250  float xm = fts.momentum().x();
251  float ym = fts.momentum().y();
252 
253  return ((x0 * xm + y0 * ym) > 0);
254 }
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
T perp() const
Definition: PV3DBase.h:69
tuple cont
load Luminosity info ##
Definition: generateEDF.py:628
std::vector< const DetLayer * > compatibleEndcapLayers(const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
static const double slope[3]
void inOutBackward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
T y() const
Definition: PV3DBase.h:60
void outInForward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
void inOutForward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
PropagationDirection
DirectMuonNavigation(const edm::ESHandle< MuonDetLayerGeometry > &)
bool outward(const FreeTrajectoryState &fts) const
std::vector< const DetLayer * > compatibleLayers(const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
bool checkCompatible(const FreeTrajectoryState &fts, const BarrelDetLayer *) const
void outInBarrel(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
T z() const
Definition: PV3DBase.h:61
void outInBackward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
GlobalVector momentum() const
GlobalPoint position() const
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.
void inOutBarrel(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
virtual const BoundDisk & specificSurface() const final
def rm
Definition: eostools.py:363
T x() const
Definition: PV3DBase.h:59