CMS 3D CMS Logo

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