CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/RecoMuon/Navigation/src/MuonNavigationSchool.cc

Go to the documentation of this file.
00001 
00020 #include "RecoMuon/Navigation/interface/MuonNavigationSchool.h"
00021 
00022 /* Collaborating Class Header */
00023 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
00024 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
00025 #include "TrackingTools/DetLayers/interface/NavigationSetter.h"
00026 #include "DataFormats/GeometrySurface/interface/BoundCylinder.h"
00027 #include "DataFormats/GeometrySurface/interface/BoundDisk.h"
00028 #include "RecoMuon/DetLayers/interface/MuonDetLayerGeometry.h"
00029 #include "RecoMuon/Navigation/interface/MuonBarrelNavigableLayer.h"
00030 #include "RecoMuon/Navigation/interface/MuonForwardNavigableLayer.h"
00031 #include "RecoMuon/Navigation/interface/MuonEtaRange.h"
00032 #include "RecoMuon/Navigation/interface/MuonDetLayerMap.h"
00033 #include "Utilities/General/interface/CMSexception.h"
00034 
00035 #include <algorithm>
00036 #include <iostream>
00037 using namespace std;
00038 
00040 MuonNavigationSchool::MuonNavigationSchool(const MuonDetLayerGeometry * muonLayout, bool enableRPC ) : theMuonDetLayerGeometry(muonLayout) {
00041 
00042   theAllDetLayersInSystem=&muonLayout->allLayers(); 
00043 
00044   // get all barrel DetLayers (DT + optional RPC) 
00045   vector<DetLayer*> barrel;
00046   if ( enableRPC ) barrel = muonLayout->allBarrelLayers();
00047   else barrel = muonLayout->allDTLayers();
00048 
00049   for ( vector<DetLayer*>::const_iterator i = barrel.begin(); i != barrel.end(); i++ ) {
00050     BarrelDetLayer* mbp = dynamic_cast<BarrelDetLayer*>(*i);
00051     if ( mbp == 0 ) throw Genexception("Bad BarrelDetLayer");
00052     addBarrelLayer(mbp);
00053   }
00054 
00055   // get all endcap DetLayers (CSC + optional RPC)
00056   vector<DetLayer*> endcap;
00057   if ( enableRPC ) endcap = muonLayout->allEndcapLayers();
00058   else endcap = muonLayout->allCSCLayers();
00059 
00060   for ( vector<DetLayer*>::const_iterator i = endcap.begin(); i != endcap.end(); i++ ) {
00061     ForwardDetLayer* mep = dynamic_cast<ForwardDetLayer*>(*i);
00062     if ( mep == 0 ) throw Genexception("Bad ForwardDetLayer");
00063     addEndcapLayer(mep);
00064   }
00065 
00066   // create outward links for all DetLayers
00067   linkBarrelLayers();
00068   linkEndcapLayers(theForwardLayers,theForwardNLC);
00069   linkEndcapLayers(theBackwardLayers,theBackwardNLC);
00070 
00071   // create inverse links
00072   createInverseLinks();
00073 
00074 }
00075 
00076 
00078 MuonNavigationSchool::~MuonNavigationSchool() {
00079 
00080    for_each(theBarrelNLC.begin(),theBarrelNLC.end(), delete_layer());
00081    for_each(theForwardNLC.begin(),theForwardNLC.end(), delete_layer());
00082    for_each(theBackwardNLC.begin(),theBackwardNLC.end(), delete_layer());
00083 
00084 }
00085 
00086 
00088 MuonNavigationSchool::StateType 
00089 MuonNavigationSchool::navigableLayers() const {
00090 
00091   StateType result;
00092   
00093   vector<MuonBarrelNavigableLayer*>::const_iterator ib;
00094   vector<MuonForwardNavigableLayer*>::const_iterator ie;
00095 
00096   for ( ib = theBarrelNLC.begin(); ib != theBarrelNLC.end(); ib++ ) {
00097     result.push_back(*ib);
00098   }
00099 
00100   for ( ie = theForwardNLC.begin(); ie != theForwardNLC.end(); ie++ ) {
00101     result.push_back(*ie);
00102   }
00103 
00104   for ( ie = theBackwardNLC.begin(); ie != theBackwardNLC.end(); ie++ ) {
00105     result.push_back(*ie);
00106   }
00107   
00108   return result;
00109 
00110 }
00111 
00112 
00114 void MuonNavigationSchool::addBarrelLayer(BarrelDetLayer* mbp) {
00115 
00116   const BoundCylinder& bc = mbp->specificSurface();
00117   float radius = bc.radius();
00118   float length = bc.bounds().length()/2.;
00119 
00120   float eta_max = calculateEta(radius, length);
00121   float eta_min = -eta_max;
00122 
00123   theBarrelLayers[mbp] = MuonEtaRange(eta_max, eta_min);
00124 
00125 }
00126 
00127 
00129 void MuonNavigationSchool::addEndcapLayer(ForwardDetLayer* mep) {
00130 
00131   const BoundDisk& bd = mep->specificSurface();
00132   float outRadius = bd.outerRadius();
00133   float inRadius = bd.innerRadius();
00134   float thick = bd.bounds().length()/2.;
00135   float z = bd.position().z();
00136 
00137   if ( z > 0. ) {
00138     float eta_min = calculateEta(outRadius, z-thick);
00139     float eta_max = calculateEta(inRadius, z+thick);
00140     theForwardLayers[mep] = MuonEtaRange(eta_max, eta_min);
00141   } else {
00142     float eta_max = calculateEta(outRadius, z+thick);
00143     float eta_min = calculateEta(inRadius, z-thick);
00144     theBackwardLayers[mep] = MuonEtaRange(eta_max, eta_min);
00145   }
00146 
00147 }
00148 
00149 
00151 float MuonNavigationSchool::calculateEta(const float& r, const float& z) const {
00152 
00153   if ( z > 0 ) return -log((tan(atan(r/z)/2.)));
00154   return log(-(tan(atan(r/z)/2.)));
00155 
00156 }
00157 
00159 void MuonNavigationSchool::linkBarrelLayers() {
00160 
00161   for (MapBI bl  = theBarrelLayers.begin();
00162              bl != theBarrelLayers.end(); bl++) {
00163 
00164     MuonEtaRange range = (*bl).second;
00165 
00166     // first add next barrel layer
00167     MapBI plusOne(bl);
00168     plusOne++;
00169     MapB outerBarrel;
00170     MapB allOuterBarrel;
00171     if ( plusOne != theBarrelLayers.end() ) { outerBarrel.insert(*plusOne);}
00172     // add all outer barrel layers
00173     for ( MapBI iMBI = plusOne; iMBI!= theBarrelLayers.end(); iMBI++){
00174       allOuterBarrel.insert(*iMBI);
00175     }
00176     // then add all compatible backward layers with an eta criteria
00177     MapE allOuterBackward;
00178     for (MapEI el  = theBackwardLayers.begin();
00179                el != theBackwardLayers.end(); el++) {
00180       if ( (*el).second.isCompatible(range) ) {
00181         allOuterBackward.insert(*el);
00182       }
00183     }
00184     //add the backward next layer with an eta criteria
00185     MapE outerBackward;
00186     for (MapEI el  = theBackwardLayers.begin();
00187                el != theBackwardLayers.end(); el++) {
00188       if ( (*el).second.isCompatible(range) ) {
00189         outerBackward.insert(*el);
00190         break;
00191       }
00192     }
00193 
00194     // then add all compatible forward layers with an eta criteria
00195     MapE allOuterForward;
00196     for (MapEI el  = theForwardLayers.begin();
00197                el != theForwardLayers.end(); el++) {
00198       if ( (*el).second.isCompatible(range) ) {
00199         allOuterForward.insert(*el);
00200       }
00201     }
00202 
00203     // then add forward next layer with an eta criteria
00204     MapE outerForward;
00205     for (MapEI el  = theForwardLayers.begin();
00206                el != theForwardLayers.end(); el++) {
00207       if ( (*el).second.isCompatible(range) ) {
00208         outerForward.insert(*el);
00209         break;
00210       }
00211     }
00212 
00213     theBarrelNLC.push_back(new MuonBarrelNavigableLayer(
00214                        (*bl).first,outerBarrel, outerBackward, outerForward,
00215                        allOuterBarrel,allOuterBackward,allOuterForward));
00216 
00217   }
00218 
00219 }
00221 void MuonNavigationSchool::linkEndcapLayers(const MapE& layers,
00222                                             vector<MuonForwardNavigableLayer*>& result) {
00223 
00224   for (MapEI el = layers.begin(); el != layers.end(); el++) {
00225 
00226     MuonEtaRange range = (*el).second;
00227     // first add next endcap layer (if compatible)
00228     MapEI plusOne(el); 
00229     plusOne++;
00230     MapE outerLayers;
00231     if ( plusOne != layers.end() && (*plusOne).second.isCompatible(range) ) {
00232         outerLayers.insert(*plusOne);
00233       if ( !range.isInside((*plusOne).second) ) {
00234         // then look if the next layer has a wider eta range, if so add it
00235         MapEI tmpel(plusOne);
00236         tmpel++;
00237         MuonEtaRange max((*plusOne).second);
00238         for ( MapEI l = tmpel; l != layers.end(); l++ ) {
00239           MuonEtaRange next = (*l).second;
00240           if ( next.isCompatible(max) && !range.isInside(next) &&
00241                !next.isInside(max) && next.subtract(max).isInside(range) ) {
00242             max = max.add(next);
00243             outerLayers.insert(*l);
00244           }
00245         }
00246       }
00247     }
00248 
00249     MapE allOuterLayers;
00250     for (MapEI iMEI = plusOne; iMEI!=layers.end(); iMEI++){
00251       if ((*iMEI).second.isCompatible(range)) allOuterLayers.insert(*iMEI);
00252     }
00253     
00254     result.push_back(new MuonForwardNavigableLayer(
00255                    (*el).first,outerLayers, allOuterLayers));
00256   }
00257 
00258 }
00259 
00260 
00262 void MuonNavigationSchool::createInverseLinks() const {
00263 
00264   // set outward link
00265   NavigationSetter setter(*this);
00266 
00267   // find for each layer which are the layers pointing to it
00268   typedef map<const DetLayer*, MapB, less<const DetLayer*> > BarrelMapType;
00269   typedef map<const DetLayer*, MapE, less<const DetLayer*> > ForwardMapType;
00270 
00271   // map of all DetLayers which can reach a specific DetLayer
00272   BarrelMapType reachedBarrelLayersMap;
00273   ForwardMapType reachedForwardLayersMap;
00274 
00275   // map of all DetLayers which is compatible with a specific DetLayer
00276   BarrelMapType compatibleBarrelLayersMap;
00277   ForwardMapType compatibleForwardLayersMap;
00278 
00279   // collect all reacheable layers starting from a barrel layer
00280   for ( MapBI bli  = theBarrelLayers.begin(); 
00281               bli != theBarrelLayers.end(); bli++ ) {
00282     // barrel
00283     MuonBarrelNavigableLayer* mbnl =
00284       dynamic_cast<MuonBarrelNavigableLayer*>(((*bli).first)->navigableLayer());
00285     MapB reacheableB = mbnl->getOuterBarrelLayers();
00286     for (MapBI i = reacheableB.begin(); i != reacheableB.end(); i++ ) {
00287       reachedBarrelLayersMap[(*i).first].insert(*bli);
00288     }
00289     MapB compatibleB = mbnl->getAllOuterBarrelLayers();
00290     for (MapBI i = compatibleB.begin(); i != compatibleB.end(); i++ ) {
00291       compatibleBarrelLayersMap[(*i).first].insert(*bli);
00292     }
00293     MapE reacheableE = mbnl->getOuterBackwardLayers();
00294     for (MapEI i = reacheableE.begin(); i != reacheableE.end(); i++ ) {
00295       reachedBarrelLayersMap[(*i).first].insert(*bli);
00296     }
00297     reacheableE = mbnl->getOuterForwardLayers();
00298     for (MapEI i = reacheableE.begin(); i != reacheableE.end(); i++ ) {
00299       reachedBarrelLayersMap[(*i).first].insert(*bli);
00300     }
00301     MapE compatibleE = mbnl->getAllOuterBackwardLayers();
00302     for (MapEI i = compatibleE.begin(); i != compatibleE.end(); i++ ) {
00303       compatibleBarrelLayersMap[(*i).first].insert(*bli);
00304     }
00305     compatibleE = mbnl->getAllOuterForwardLayers();
00306     for (MapEI i = compatibleE.begin(); i != compatibleE.end(); i++ ) {
00307       compatibleBarrelLayersMap[(*i).first].insert(*bli);
00308     }
00309 
00310   }
00311 
00312   // collect all reacheable layer starting from a backward layer
00313   for ( MapEI eli  = theBackwardLayers.begin(); 
00314               eli != theBackwardLayers.end(); eli++ ) {
00315     MapE reacheableE =
00316       dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer())->getOuterEndcapLayers();
00317     for (MapEI i = reacheableE.begin(); i != reacheableE.end(); i++ ) {
00318       reachedForwardLayersMap[(*i).first].insert(*eli);
00319     }
00320   // collect all compatible layer starting from a backward layer
00321     MapE compatibleE =
00322       dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer())->getAllOuterEndcapLayers();
00323     for (MapEI i = compatibleE.begin(); i != compatibleE.end(); i++ ) {
00324       compatibleForwardLayersMap[(*i).first].insert(*eli);
00325     }
00326   }
00327 
00328   for ( MapEI eli  = theForwardLayers.begin(); 
00329               eli != theForwardLayers.end(); eli++ ) {
00330   // collect all reacheable layer starting from a forward layer
00331     MapE reacheableE =
00332       dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer())->getOuterEndcapLayers();
00333     for (MapEI i = reacheableE.begin(); i != reacheableE.end(); i++ ) {
00334       reachedForwardLayersMap[(*i).first].insert(*eli);
00335     }
00336   // collect all compatible layer starting from a forward layer
00337     MapE compatibleE =
00338       dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer())->getAllOuterEndcapLayers();
00339     for (MapEI i = compatibleE.begin(); i != compatibleE.end(); i++ ) {
00340       compatibleForwardLayersMap[(*i).first].insert(*eli);
00341     }
00342   }
00343 
00344   // now set inverse link for barrel layers
00345   for ( MapBI bli  = theBarrelLayers.begin(); 
00346               bli != theBarrelLayers.end(); bli++ ) {
00347     MuonBarrelNavigableLayer* mbnl =
00348       dynamic_cast<MuonBarrelNavigableLayer*>(((*bli).first)->navigableLayer());
00349     mbnl->setInwardLinks(reachedBarrelLayersMap[(*bli).first]);
00350     mbnl->setInwardCompatibleLinks(compatibleBarrelLayersMap[(*bli).first]);
00351 
00352   }
00353   //BACKWARD
00354   for ( MapEI eli  = theBackwardLayers.begin(); 
00355               eli != theBackwardLayers.end(); eli++ ) {
00356     MuonForwardNavigableLayer* mfnl =      
00357       dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer());
00358     // for backward next layers
00359     mfnl->setInwardLinks(reachedBarrelLayersMap[(*eli).first],
00360                          reachedForwardLayersMap[(*eli).first]);
00361   // for backward compatible layers
00362     mfnl->setInwardCompatibleLinks(compatibleBarrelLayersMap[(*eli).first],
00363                          compatibleForwardLayersMap[(*eli).first]);
00364   }
00365   //FORWARD
00366   for ( MapEI eli  = theForwardLayers.begin(); 
00367               eli != theForwardLayers.end(); eli++ ) {
00368     MuonForwardNavigableLayer* mfnl = 
00369       dynamic_cast<MuonForwardNavigableLayer*>(((*eli).first)->navigableLayer());
00370   // and for forward next layers
00371     mfnl->setInwardLinks(reachedBarrelLayersMap[(*eli).first],
00372                          reachedForwardLayersMap[(*eli).first]);
00373   // and for forward compatible layers
00374     mfnl->setInwardCompatibleLinks(compatibleBarrelLayersMap[(*eli).first],
00375                          compatibleForwardLayersMap[(*eli).first]);
00376   }
00377 
00378 }