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