CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonNavigationSchool.cc
Go to the documentation of this file.
1 
19 
20 /* Collaborating Class Header */
32 
33 #include <algorithm>
34 #include <iostream>
35 using namespace std;
36 
38 MuonNavigationSchool::MuonNavigationSchool(const MuonDetLayerGeometry * muonLayout, bool enableRPC ) : theMuonDetLayerGeometry(muonLayout) {
39 
40  theAllDetLayersInSystem=&muonLayout->allLayers();
41 
42  // get all barrel DetLayers (DT + optional RPC)
43  vector<DetLayer*> barrel;
44  if ( enableRPC ) barrel = muonLayout->allBarrelLayers();
45  else barrel = muonLayout->allDTLayers();
46 
47  for ( vector<DetLayer*>::const_iterator i = barrel.begin(); i != barrel.end(); i++ ) {
48  BarrelDetLayer* mbp = dynamic_cast<BarrelDetLayer*>(*i);
49  if ( mbp == 0 ) throw cms::Exception("MuonNavigationSchool", "Bad BarrelDetLayer");
50  addBarrelLayer(mbp);
51  }
52 
53  // get all endcap DetLayers (CSC + optional RPC)
54  vector<DetLayer*> endcap;
55  if ( enableRPC ) endcap = muonLayout->allEndcapLayers();
56  else endcap = muonLayout->allCSCLayers();
57 
58  for ( vector<DetLayer*>::const_iterator i = endcap.begin(); i != endcap.end(); i++ ) {
59  ForwardDetLayer* mep = dynamic_cast<ForwardDetLayer*>(*i);
60  if ( mep == 0 ) throw cms::Exception("MuonNavigationSchool", "Bad ForwardDetLayer");
61  addEndcapLayer(mep);
62  }
63 
64  // create outward links for all DetLayers
68 
69  // create inverse links
71 
72 }
73 
74 
77 
78  for_each(theBarrelNLC.begin(),theBarrelNLC.end(), delete_layer());
79  for_each(theForwardNLC.begin(),theForwardNLC.end(), delete_layer());
80  for_each(theBackwardNLC.begin(),theBackwardNLC.end(), delete_layer());
81 
82 }
83 
84 
88 
90 
91  vector<MuonBarrelNavigableLayer*>::const_iterator ib;
92  vector<MuonForwardNavigableLayer*>::const_iterator ie;
93 
94  for ( ib = theBarrelNLC.begin(); ib != theBarrelNLC.end(); ib++ ) {
95  result.push_back(*ib);
96  }
97 
98  for ( ie = theForwardNLC.begin(); ie != theForwardNLC.end(); ie++ ) {
99  result.push_back(*ie);
100  }
101 
102  for ( ie = theBackwardNLC.begin(); ie != theBackwardNLC.end(); ie++ ) {
103  result.push_back(*ie);
104  }
105 
106  return result;
107 
108 }
109 
110 
113 
114  const BoundCylinder& bc = mbp->specificSurface();
115  float radius = bc.radius();
116  float length = bc.bounds().length()/2.;
117 
118  float eta_max = calculateEta(radius, length);
119  float eta_min = -eta_max;
120 
121  theBarrelLayers[mbp] = MuonEtaRange(eta_max, eta_min);
122 
123 }
124 
125 
128 
129  const BoundDisk& bd = mep->specificSurface();
130  float outRadius = bd.outerRadius();
131  float inRadius = bd.innerRadius();
132  float thick = bd.bounds().length()/2.;
133  float z = bd.position().z();
134 
135  if ( z > 0. ) {
136  float eta_min = calculateEta(outRadius, z-thick);
137  float eta_max = calculateEta(inRadius, z+thick);
138  theForwardLayers[mep] = MuonEtaRange(eta_max, eta_min);
139  } else {
140  float eta_max = calculateEta(outRadius, z+thick);
141  float eta_min = calculateEta(inRadius, z-thick);
142  theBackwardLayers[mep] = MuonEtaRange(eta_max, eta_min);
143  }
144 
145 }
146 
147 
149 float MuonNavigationSchool::calculateEta(const float& r, const float& z) const {
150 
151  if ( z > 0 ) return -log((tan(atan(r/z)/2.)));
152  return log(-(tan(atan(r/z)/2.)));
153 
154 }
155 
158 
159  for (MapBI bl = theBarrelLayers.begin();
160  bl != theBarrelLayers.end(); bl++) {
161 
162  MuonEtaRange range = (*bl).second;
163 
164  // first add next barrel layer
165  MapBI plusOne(bl);
166  plusOne++;
167  MapB outerBarrel;
168  MapB allOuterBarrel;
169  if ( plusOne != theBarrelLayers.end() ) { outerBarrel.insert(*plusOne);}
170  // add all outer barrel layers
171  for ( MapBI iMBI = plusOne; iMBI!= theBarrelLayers.end(); iMBI++){
172  allOuterBarrel.insert(*iMBI);
173  }
174  // then add all compatible backward layers with an eta criteria
175  MapE allOuterBackward;
176  for (MapEI el = theBackwardLayers.begin();
177  el != theBackwardLayers.end(); el++) {
178  if ( (*el).second.isCompatible(range) ) {
179  allOuterBackward.insert(*el);
180  }
181  }
182  //add the backward next layer with an eta criteria
183  MapE outerBackward;
184  for (MapEI el = theBackwardLayers.begin();
185  el != theBackwardLayers.end(); el++) {
186  if ( (*el).second.isCompatible(range) ) {
187  outerBackward.insert(*el);
188  break;
189  }
190  }
191 
192  // then add all compatible forward layers with an eta criteria
193  MapE allOuterForward;
194  for (MapEI el = theForwardLayers.begin();
195  el != theForwardLayers.end(); el++) {
196  if ( (*el).second.isCompatible(range) ) {
197  allOuterForward.insert(*el);
198  }
199  }
200 
201  // then add forward next layer with an eta criteria
202  MapE outerForward;
203  for (MapEI el = theForwardLayers.begin();
204  el != theForwardLayers.end(); el++) {
205  if ( (*el).second.isCompatible(range) ) {
206  outerForward.insert(*el);
207  break;
208  }
209  }
210 
212  (*bl).first,outerBarrel, outerBackward, outerForward,
213  allOuterBarrel,allOuterBackward,allOuterForward));
214 
215  }
216 
217 }
220  vector<MuonForwardNavigableLayer*>& result) {
221 
222  for (MapEI el = layers.begin(); el != layers.end(); el++) {
223 
224  MuonEtaRange range = (*el).second;
225  // first add next endcap layer (if compatible)
226  MapEI plusOne(el);
227  plusOne++;
228  MapE outerLayers;
229  if ( plusOne != layers.end() && (*plusOne).second.isCompatible(range) ) {
230  outerLayers.insert(*plusOne);
231  if ( !range.isInside((*plusOne).second) ) {
232  // then look if the next layer has a wider eta range, if so add it
233  MapEI tmpel(plusOne);
234  tmpel++;
235  MuonEtaRange max((*plusOne).second);
236  for ( MapEI l = tmpel; l != layers.end(); l++ ) {
237  MuonEtaRange next = (*l).second;
238  if ( next.isCompatible(max) && !range.isInside(next) &&
239  !next.isInside(max) && next.subtract(max).isInside(range) ) {
240  max = max.add(next);
241  outerLayers.insert(*l);
242  }
243  }
244  }
245  }
246 
247  MapE allOuterLayers;
248  for (MapEI iMEI = plusOne; iMEI!=layers.end(); iMEI++){
249  if ((*iMEI).second.isCompatible(range)) allOuterLayers.insert(*iMEI);
250  }
251 
252  result.push_back(new MuonForwardNavigableLayer(
253  (*el).first,outerLayers, allOuterLayers));
254  }
255 
256 }
257 
258 
261 
262  // set outward link
263  NavigationSetter setter(*this);
264 
265  // find for each layer which are the layers pointing to it
266  typedef map<const DetLayer*, MapB, less<const DetLayer*> > BarrelMapType;
267  typedef map<const DetLayer*, MapE, less<const DetLayer*> > ForwardMapType;
268 
269  // map of all DetLayers which can reach a specific DetLayer
270  BarrelMapType reachedBarrelLayersMap;
271  ForwardMapType reachedForwardLayersMap;
272 
273  // map of all DetLayers which is compatible with a specific DetLayer
274  BarrelMapType compatibleBarrelLayersMap;
275  ForwardMapType compatibleForwardLayersMap;
276 
277  // collect all reacheable layers starting from a barrel layer
278  for ( MapBI bli = theBarrelLayers.begin();
279  bli != theBarrelLayers.end(); bli++ ) {
280  // barrel
282  dynamic_cast<MuonBarrelNavigableLayer*>(((*bli).first)->navigableLayer());
283  MapB reacheableB = mbnl->getOuterBarrelLayers();
284  for (MapBI i = reacheableB.begin(); i != reacheableB.end(); i++ ) {
285  reachedBarrelLayersMap[(*i).first].insert(*bli);
286  }
287  MapB compatibleB = mbnl->getAllOuterBarrelLayers();
288  for (MapBI i = compatibleB.begin(); i != compatibleB.end(); i++ ) {
289  compatibleBarrelLayersMap[(*i).first].insert(*bli);
290  }
291  MapE reacheableE = mbnl->getOuterBackwardLayers();
292  for (MapEI i = reacheableE.begin(); i != reacheableE.end(); i++ ) {
293  reachedBarrelLayersMap[(*i).first].insert(*bli);
294  }
295  reacheableE = mbnl->getOuterForwardLayers();
296  for (MapEI i = reacheableE.begin(); i != reacheableE.end(); i++ ) {
297  reachedBarrelLayersMap[(*i).first].insert(*bli);
298  }
299  MapE compatibleE = mbnl->getAllOuterBackwardLayers();
300  for (MapEI i = compatibleE.begin(); i != compatibleE.end(); i++ ) {
301  compatibleBarrelLayersMap[(*i).first].insert(*bli);
302  }
303  compatibleE = mbnl->getAllOuterForwardLayers();
304  for (MapEI i = compatibleE.begin(); i != compatibleE.end(); i++ ) {
305  compatibleBarrelLayersMap[(*i).first].insert(*bli);
306  }
307 
308  }
309 
310  // collect all reacheable layer starting from a backward layer
311  for ( MapEI eli = theBackwardLayers.begin();
312  eli != theBackwardLayers.end(); eli++ ) {
313  MapE reacheableE =
314  dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer())->getOuterEndcapLayers();
315  for (MapEI i = reacheableE.begin(); i != reacheableE.end(); i++ ) {
316  reachedForwardLayersMap[(*i).first].insert(*eli);
317  }
318  // collect all compatible layer starting from a backward layer
319  MapE compatibleE =
320  dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer())->getAllOuterEndcapLayers();
321  for (MapEI i = compatibleE.begin(); i != compatibleE.end(); i++ ) {
322  compatibleForwardLayersMap[(*i).first].insert(*eli);
323  }
324  }
325 
326  for ( MapEI eli = theForwardLayers.begin();
327  eli != theForwardLayers.end(); eli++ ) {
328  // collect all reacheable layer starting from a forward layer
329  MapE reacheableE =
330  dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer())->getOuterEndcapLayers();
331  for (MapEI i = reacheableE.begin(); i != reacheableE.end(); i++ ) {
332  reachedForwardLayersMap[(*i).first].insert(*eli);
333  }
334  // collect all compatible layer starting from a forward layer
335  MapE compatibleE =
336  dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer())->getAllOuterEndcapLayers();
337  for (MapEI i = compatibleE.begin(); i != compatibleE.end(); i++ ) {
338  compatibleForwardLayersMap[(*i).first].insert(*eli);
339  }
340  }
341 
342  // now set inverse link for barrel layers
343  for ( MapBI bli = theBarrelLayers.begin();
344  bli != theBarrelLayers.end(); bli++ ) {
346  dynamic_cast<MuonBarrelNavigableLayer*>(((*bli).first)->navigableLayer());
347  mbnl->setInwardLinks(reachedBarrelLayersMap[(*bli).first]);
348  mbnl->setInwardCompatibleLinks(compatibleBarrelLayersMap[(*bli).first]);
349 
350  }
351  //BACKWARD
352  for ( MapEI eli = theBackwardLayers.begin();
353  eli != theBackwardLayers.end(); eli++ ) {
355  dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer());
356  // for backward next layers
357  mfnl->setInwardLinks(reachedBarrelLayersMap[(*eli).first],
358  reachedForwardLayersMap[(*eli).first]);
359  // for backward compatible layers
360  mfnl->setInwardCompatibleLinks(compatibleBarrelLayersMap[(*eli).first],
361  compatibleForwardLayersMap[(*eli).first]);
362  }
363  //FORWARD
364  for ( MapEI eli = theForwardLayers.begin();
365  eli != theForwardLayers.end(); eli++ ) {
367  dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer());
368  // and for forward next layers
369  mfnl->setInwardLinks(reachedBarrelLayersMap[(*eli).first],
370  reachedForwardLayersMap[(*eli).first]);
371  // and for forward compatible layers
372  mfnl->setInwardCompatibleLinks(compatibleBarrelLayersMap[(*eli).first],
373  compatibleForwardLayersMap[(*eli).first]);
374  }
375 
376 }
std::vector< MuonForwardNavigableLayer * > theForwardNLC
int i
Definition: DBlmapReader.cc:9
MuonEtaRange add(const MuonEtaRange &) const
create maximum of ranges
Definition: MuonEtaRange.cc:59
int ib
Definition: cuy.py:660
MapE theBackwardLayers
+z endcap
const std::vector< DetLayer * > & allCSCLayers() const
return the CSC DetLayers (endcap), -Z to +Z
std::map< ForwardDetLayer *, MuonEtaRange, MuonDetLayerComp > MapE
MapB::const_iterator MapBI
const std::vector< DetLayer * > & allDTLayers() const
return the DT DetLayers (barrel), inside-out
float float float z
MuonNavigationSchool(const MuonDetLayerGeometry *, bool enableRPC=true)
Constructor.
void setInwardLinks(const MapB &, const MapE &)
set inward links
std::map< BarrelDetLayer *, MuonEtaRange, MuonDetLayerComp > MapB
void createInverseLinks() const
establish inward links
void setInwardCompatibleLinks(const MapB &)
std::vector< MuonForwardNavigableLayer * > theBackwardNLC
float calculateEta(const float &, const float &) const
calculate pseudorapidity from r and z
MapE::const_iterator MapEI
const T & max(const T &a, const T &b)
std::vector< MuonBarrelNavigableLayer * > theBarrelNLC
-z endcap
tuple result
Definition: query.py:137
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
void linkBarrelLayers()
link barrel layers
std::vector< NavigableLayer * > StateType
const std::vector< DetLayer * > & allEndcapLayers() const
return all endcap DetLayers (CSC+RPC), -Z to +Z
bool isInside(float eta, float error=0.) const
Definition: MuonEtaRange.cc:43
const std::vector< DetLayer * > & allLayers() const
return all layers (DT+CSC+RPC), order: backward, barrel, forward
virtual const BoundDisk & specificSurface() const GCC11_FINAL
void linkEndcapLayers(const MapE &, std::vector< MuonForwardNavigableLayer * > &)
link endcap layers
~MuonNavigationSchool()
Destructor.
bool isCompatible(const MuonEtaRange &range) const
true if this overlaps with range
Definition: MuonEtaRange.cc:54
void setInwardLinks(const MapB &)
set inward links
void setInwardCompatibleLinks(const MapB &, const MapE &)
void addEndcapLayer(ForwardDetLayer *)
add endcap layer (backward and forward)
MuonEtaRange subtract(const MuonEtaRange &) const
create new range of size this minus range
Definition: MuonEtaRange.cc:65
const std::vector< DetLayer * > * theAllDetLayersInSystem
virtual const BoundCylinder & specificSurface() const GCC11_FINAL
Extension of the interface.
void addBarrelLayer(BarrelDetLayer *)
add barrel layer
const std::vector< DetLayer * > & allBarrelLayers() const
return all barrel DetLayers (DT+RPC), inside-out
virtual StateType navigableLayers() const
return navigable layers, from base class