CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
DirectMuonNavigation Class Reference

#include <DirectMuonNavigation.h>

Public Member Functions

DirectMuonNavigationclone () const
 
std::vector< const DetLayer * > compatibleEndcapLayers (const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
 
std::vector< const DetLayer * > compatibleLayers (const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
 
 DirectMuonNavigation (const edm::ESHandle< MuonDetLayerGeometry > &)
 
 DirectMuonNavigation (const edm::ESHandle< MuonDetLayerGeometry > &, const edm::ParameterSet &)
 
 ~DirectMuonNavigation ()
 

Private Member Functions

bool checkCompatible (const FreeTrajectoryState &fts, const BarrelDetLayer *) const
 
bool checkCompatible (const FreeTrajectoryState &fts, const ForwardDetLayer *) const
 
void inOutBackward (const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
 
void inOutBarrel (const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
 
void inOutForward (const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
 
void outInBackward (const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
 
void outInBarrel (const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
 
void outInForward (const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
 
bool outward (const FreeTrajectoryState &fts) const
 

Private Attributes

float epsilon_
 
bool theBarrelFlag
 
bool theEndcapFlag
 
edm::ESHandle
< MuonDetLayerGeometry
theMuonDetLayerGeometry
 

Detailed Description

Definition at line 22 of file DirectMuonNavigation.h.

Constructor & Destructor Documentation

DirectMuonNavigation::DirectMuonNavigation ( const edm::ESHandle< MuonDetLayerGeometry > &  muonLayout)

Definition at line 21 of file DirectMuonNavigation.cc.

Referenced by clone().

21  : theMuonDetLayerGeometry(muonLayout), epsilon_(100.), theEndcapFlag(true), theBarrelFlag(true) {
22 }
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
DirectMuonNavigation::DirectMuonNavigation ( const edm::ESHandle< MuonDetLayerGeometry > &  muonLayout,
const edm::ParameterSet par 
)

Definition at line 24 of file DirectMuonNavigation.cc.

24  : theMuonDetLayerGeometry(muonLayout), epsilon_(100.), theEndcapFlag(par.getParameter<bool>("Endcap")), theBarrelFlag(par.getParameter<bool>("Barrel")) {
25 
26 }
T getParameter(std::string const &) const
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
DirectMuonNavigation::~DirectMuonNavigation ( )
inline

Definition at line 36 of file DirectMuonNavigation.h.

36 {}

Member Function Documentation

bool DirectMuonNavigation::checkCompatible ( const FreeTrajectoryState fts,
const BarrelDetLayer dl 
) const
private

Definition at line 206 of file DirectMuonNavigation.cc.

References BoundSurface::bounds(), epsilon_, Bounds::length(), FreeTrajectoryState::momentum(), outward(), PV3DBase< T, PVType, FrameType >::perp(), FreeTrajectoryState::position(), Cylinder::radius(), CosmicsPD_Skims::radius, submit::rm, slope, BarrelDetLayer::specificSurface(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by inOutBackward(), inOutBarrel(), inOutForward(), outInBackward(), outInBarrel(), and outInForward().

206  {
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 }
virtual float length() const =0
T perp() const
Definition: PV3DBase.h:66
static const double slope[3]
bool outward(const FreeTrajectoryState &fts) const
Scalar radius() const
Radius of the cylinder.
Definition: Cylinder.h:55
T z() const
Definition: PV3DBase.h:58
string rm
Definition: submit.py:76
GlobalVector momentum() const
GlobalPoint position() const
const Bounds & bounds() const
Definition: BoundSurface.h:89
virtual const BoundCylinder & specificSurface() const
Extension of the interface.
bool DirectMuonNavigation::checkCompatible ( const FreeTrajectoryState fts,
const ForwardDetLayer dl 
) const
private

Definition at line 223 of file DirectMuonNavigation.cc.

References epsilon_, BoundDisk::innerRadius(), FreeTrajectoryState::momentum(), BoundDisk::outerRadius(), outward(), PV3DBase< T, PVType, FrameType >::perp(), GloballyPositioned< T >::position(), FreeTrajectoryState::position(), submit::rm, slope, ForwardDetLayer::specificSurface(), detailsBasic3DVector::z, and PV3DBase< T, PVType, FrameType >::z().

223  {
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 }
T perp() const
Definition: PV3DBase.h:66
static const double slope[3]
double double double z
bool outward(const FreeTrajectoryState &fts) const
T z() const
Definition: PV3DBase.h:58
virtual const BoundDisk & specificSurface() const
string rm
Definition: submit.py:76
GlobalVector momentum() const
GlobalPoint position() const
float outerRadius() const
The outer radius of the disk.
Definition: BoundDisk.cc:10
float innerRadius() const
The inner radius of the disk.
Definition: BoundDisk.cc:6
const PositionType & position() const
DirectMuonNavigation* DirectMuonNavigation::clone ( void  ) const
inline

Definition at line 31 of file DirectMuonNavigation.h.

References DirectMuonNavigation().

31  {
32  return new DirectMuonNavigation(*this);
33  }
DirectMuonNavigation(const edm::ESHandle< MuonDetLayerGeometry > &)
vector< const DetLayer * > DirectMuonNavigation::compatibleEndcapLayers ( const FreeTrajectoryState fts,
PropagationDirection  timeDirection 
) const

Definition at line 83 of file DirectMuonNavigation.cc.

References alongMomentum, inOutForward(), FreeTrajectoryState::momentum(), oppositeToMomentum, outInBackward(), convertSQLitetoXML_cfg::output, and PV3DBase< T, PVType, FrameType >::z().

Referenced by CosmicMuonTrajectoryBuilder::build(), and CosmicMuonTrajectoryBuilder::trajectories().

84  {
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 }
void inOutForward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
T z() const
Definition: PV3DBase.h:58
void outInBackward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
GlobalVector momentum() const
dbl *** dir
Definition: mlp_gen.cc:35
vector< const DetLayer * > DirectMuonNavigation::compatibleLayers ( const FreeTrajectoryState fts,
PropagationDirection  timeDirection 
) const

Definition at line 31 of file DirectMuonNavigation.cc.

References inOutBackward(), inOutBarrel(), inOutForward(), FreeTrajectoryState::momentum(), oppositeToMomentum, outInBackward(), outInBarrel(), outInForward(), convertSQLitetoXML_cfg::output, outward(), FreeTrajectoryState::position(), theBarrelFlag, theEndcapFlag, and PV3DBase< T, PVType, FrameType >::z().

Referenced by CosmicMuonTrajectoryBuilder::build(), DynamicTruncation::compatibleDets(), DTChamberEfficiency::compatibleLayers(), StandAloneMuonFilter::compatibleLayers(), StandAloneMuonTrajectoryBuilder::propagateTheSeedTSOS(), and CosmicMuonTrajectoryBuilder::trajectories().

32  {
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 }
void inOutBackward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
void outInForward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
void inOutForward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
bool outward(const FreeTrajectoryState &fts) const
void outInBarrel(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
T z() const
Definition: PV3DBase.h:58
void outInBackward(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
GlobalVector momentum() const
GlobalPoint position() const
void inOutBarrel(const FreeTrajectoryState &, std::vector< const DetLayer * > &) const
dbl *** dir
Definition: mlp_gen.cc:35
void DirectMuonNavigation::inOutBackward ( const FreeTrajectoryState fts,
std::vector< const DetLayer * > &  output 
) const
private

Definition at line 171 of file DirectMuonNavigation.cc.

References checkCompatible(), cont, and theMuonDetLayerGeometry.

Referenced by compatibleLayers().

171  {
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 }
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
bool checkCompatible(const FreeTrajectoryState &fts, const BarrelDetLayer *) const
int cont
void DirectMuonNavigation::inOutBarrel ( const FreeTrajectoryState fts,
std::vector< const DetLayer * > &  output 
) const
private

Definition at line 102 of file DirectMuonNavigation.cc.

References Reference_intrackfit_cff::barrel, checkCompatible(), cont, and theMuonDetLayerGeometry.

Referenced by compatibleLayers().

102  {
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 }
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
bool checkCompatible(const FreeTrajectoryState &fts, const BarrelDetLayer *) const
int cont
void DirectMuonNavigation::inOutForward ( const FreeTrajectoryState fts,
std::vector< const DetLayer * > &  output 
) const
private

Definition at line 138 of file DirectMuonNavigation.cc.

References checkCompatible(), cont, and theMuonDetLayerGeometry.

Referenced by compatibleEndcapLayers(), and compatibleLayers().

138  {
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 }
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
bool checkCompatible(const FreeTrajectoryState &fts, const BarrelDetLayer *) const
int cont
void DirectMuonNavigation::outInBackward ( const FreeTrajectoryState fts,
std::vector< const DetLayer * > &  output 
) const
private

Definition at line 185 of file DirectMuonNavigation.cc.

References checkCompatible(), cont, and theMuonDetLayerGeometry.

Referenced by compatibleEndcapLayers(), and compatibleLayers().

185  {
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 }
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
bool checkCompatible(const FreeTrajectoryState &fts, const BarrelDetLayer *) const
int cont
void DirectMuonNavigation::outInBarrel ( const FreeTrajectoryState fts,
std::vector< const DetLayer * > &  output 
) const
private

Definition at line 118 of file DirectMuonNavigation.cc.

References Reference_intrackfit_cff::barrel, checkCompatible(), cont, and theMuonDetLayerGeometry.

Referenced by compatibleLayers().

118  {
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 }
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
bool checkCompatible(const FreeTrajectoryState &fts, const BarrelDetLayer *) const
int cont
void DirectMuonNavigation::outInForward ( const FreeTrajectoryState fts,
std::vector< const DetLayer * > &  output 
) const
private

Definition at line 152 of file DirectMuonNavigation.cc.

References checkCompatible(), cont, and theMuonDetLayerGeometry.

Referenced by compatibleLayers().

152  {
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 }
edm::ESHandle< MuonDetLayerGeometry > theMuonDetLayerGeometry
bool checkCompatible(const FreeTrajectoryState &fts, const BarrelDetLayer *) const
int cont
bool DirectMuonNavigation::outward ( const FreeTrajectoryState fts) const
private

Definition at line 244 of file DirectMuonNavigation.cc.

References FreeTrajectoryState::momentum(), FreeTrajectoryState::position(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by checkCompatible(), and compatibleLayers().

244  {
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 }
T y() const
Definition: PV3DBase.h:57
GlobalVector momentum() const
GlobalPoint position() const
T x() const
Definition: PV3DBase.h:56

Member Data Documentation

float DirectMuonNavigation::epsilon_
private

Definition at line 63 of file DirectMuonNavigation.h.

Referenced by checkCompatible().

bool DirectMuonNavigation::theBarrelFlag
private

Definition at line 65 of file DirectMuonNavigation.h.

Referenced by compatibleLayers().

bool DirectMuonNavigation::theEndcapFlag
private

Definition at line 64 of file DirectMuonNavigation.h.

Referenced by compatibleLayers().

edm::ESHandle<MuonDetLayerGeometry> DirectMuonNavigation::theMuonDetLayerGeometry
private