CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoMuon/GlobalTrackingTools/src/DirectTrackerNavigation.cc

Go to the documentation of this file.
00001 
00011 #include "RecoMuon/GlobalTrackingTools/interface/DirectTrackerNavigation.h"
00012 
00013 //---------------
00014 // C++ Headers --
00015 //---------------
00016 
00017 #include <algorithm>
00018 
00019 //-------------------------------
00020 // Collaborating Class Headers --
00021 //-------------------------------
00022 
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
00025 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
00026 #include "DataFormats/GeometrySurface/interface/BoundCylinder.h"
00027 #include "DataFormats/GeometrySurface/interface/BoundDisk.h"
00028 
00029 using namespace std;
00030 
00031 //
00032 // constructor
00033 //
00034 DirectTrackerNavigation::DirectTrackerNavigation(const edm::ESHandle<GeometricSearchTracker>& tkLayout, 
00035                                                  bool outOnly) : 
00036  theGeometricSearchTracker(tkLayout), theOutLayerOnlyFlag(outOnly), theEpsilon(-0.01) {
00037 
00038 }
00039 
00040 
00041 //
00042 // return compatible layers for a given trajectory state 
00043 //
00044 vector<const DetLayer*> 
00045 DirectTrackerNavigation::compatibleLayers(const FreeTrajectoryState& fts,
00046                                           PropagationDirection dir) const {
00047 
00048   bool inOut = outward(fts);
00049   double eta0 = fts.position().eta();
00050 
00051   vector<const DetLayer*> output;
00052 
00053   // check eta of DetLayers for compatibility
00054 
00055   if (inOut) { 
00056 
00057       if ( !theOutLayerOnlyFlag ) {
00058          inOutPx(fts,output);
00059 
00060          if ( eta0 > 1.55) inOutFPx(fts,output);
00061          else if ( eta0 < -1.55 ) inOutBPx(fts,output);
00062 
00063          if ( fabs(eta0) < 1.67 ) inOutTIB(fts,output); 
00064 
00065          if ( eta0 > 1.17 ) inOutFTID(fts,output);
00066          else if ( eta0 < -1.17 ) inOutBTID(fts,output);
00067       }
00068 
00069       if ( fabs(eta0) < 1.35 ) inOutTOB(fts,output);
00070 
00071       if ( eta0 > 0.97 ) inOutFTEC(fts,output);
00072       else if ( eta0 < -0.97 )  inOutBTEC(fts,output);
00073 
00074    } else {
00075       LogTrace("Muon|RecoMuon|DirectionTrackerNavigation")<<"No implementation for inward state at this moment. ";
00076    }
00077 
00078   if ( dir == oppositeToMomentum ) std::reverse(output.begin(),output.end());
00079 
00080   return output;
00081 
00082 }
00083 
00084 
00085 //
00086 //
00087 //
00088 void DirectTrackerNavigation::inOutPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00089 
00090   const vector<BarrelDetLayer*>& barrel = theGeometricSearchTracker->pixelBarrelLayers();
00091 
00092   for (vector<BarrelDetLayer*>::const_iterator iter_B = barrel.begin(); iter_B != barrel.end(); iter_B++) {
00093     if ( checkCompatible(fts,(*iter_B)) ) output.push_back((*iter_B));
00094   }
00095 
00096 }
00097 
00098 
00099 //
00100 //
00101 //
00102 void DirectTrackerNavigation::inOutTIB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00103 
00104   const vector<BarrelDetLayer*>& barrel = theGeometricSearchTracker->tibLayers();
00105 
00106   for (vector<BarrelDetLayer*>::const_iterator iter_B = barrel.begin(); iter_B != barrel.end(); iter_B++) {
00107     if ( checkCompatible(fts,(*iter_B)) ) output.push_back((*iter_B));
00108   }
00109 
00110 }
00111 
00112 
00113 //
00114 //
00115 //
00116 void DirectTrackerNavigation::inOutTOB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00117 
00118   const vector<BarrelDetLayer*>& barrel = theGeometricSearchTracker->tobLayers();
00119 
00120   for (vector<BarrelDetLayer*>::const_iterator iter_B = barrel.begin(); iter_B != barrel.end(); iter_B++) {
00121     if ( checkCompatible(fts,(*iter_B)) ) output.push_back((*iter_B));
00122   }
00123 
00124 }
00125 
00126 
00127 //
00128 //
00129 //
00130 void DirectTrackerNavigation::inOutFPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00131 
00132   const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->posPixelForwardLayers();
00133 
00134   for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
00135     if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
00136   }
00137 
00138 }
00139 
00140 
00141 //
00142 //
00143 //
00144 void DirectTrackerNavigation::inOutFTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00145 
00146   const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->posTidLayers();
00147 
00148   for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
00149     if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
00150   }
00151 
00152 }
00153 
00154 
00155 //
00156 //
00157 //
00158 void DirectTrackerNavigation::inOutFTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00159 
00160   const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->posTecLayers();
00161 
00162   for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
00163     if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
00164   }
00165 
00166 }
00167 
00168 
00169 //
00170 //
00171 //
00172 void DirectTrackerNavigation::inOutBPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00173 
00174   const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->negPixelForwardLayers();
00175 
00176   for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
00177     if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
00178   }
00179 
00180 }
00181 
00182 
00183 //
00184 //
00185 //
00186 void DirectTrackerNavigation::inOutBTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00187 
00188   const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->negTidLayers();
00189 
00190   for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
00191       if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
00192   }
00193 
00194 }
00195 
00196 
00197 //
00198 //
00199 //
00200 void DirectTrackerNavigation::inOutBTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
00201 
00202   const vector<ForwardDetLayer*>& forward = theGeometricSearchTracker->negTecLayers();
00203 
00204   for (vector<ForwardDetLayer*>::const_iterator iter_E = forward.begin(); iter_E != forward.end(); iter_E++) {
00205     if ( checkCompatible(fts,(*iter_E)) ) output.push_back((*iter_E));
00206   }
00207 
00208 }
00209 
00210 
00211 //
00212 //
00213 //
00214 bool DirectTrackerNavigation::checkCompatible(const FreeTrajectoryState& fts,const BarrelDetLayer* dl) const {
00215 
00216   float eta0 = fts.position().eta();
00217 
00218   const BoundCylinder& bc = dl->specificSurface();
00219   float radius = bc.radius();
00220   float length = bc.bounds().length()/2.;
00221 
00222   float eta = calculateEta(radius, length);
00223 
00224   return ( fabs(eta0) <= (fabs(eta) + theEpsilon) );
00225 
00226 }
00227 
00228 
00229 //
00230 //
00231 //
00232 bool DirectTrackerNavigation::checkCompatible(const FreeTrajectoryState& fts,const ForwardDetLayer* dl) const {
00233 
00234   float eta0 = fts.position().eta();
00235 
00236   const BoundDisk& bd = dl->specificSurface();
00237 
00238   float outRadius = bd.outerRadius();
00239   float inRadius = bd.innerRadius();
00240   float z = bd.position().z();
00241 
00242   float etaOut = calculateEta(outRadius, z);
00243   float etaIn = calculateEta(inRadius, z);
00244  
00245   if ( eta0 > 0 ) return ( eta0 > ( etaOut - theEpsilon) && eta0 < (etaIn + theEpsilon) );
00246   else return ( eta0 < (etaOut + theEpsilon ) && eta0 > ( etaIn - theEpsilon ) );
00247 
00248 }
00249 
00250 
00251 //
00252 //
00253 //
00254 bool DirectTrackerNavigation::outward(const FreeTrajectoryState& fts) const {
00255  
00256   return (fts.position().basicVector().dot(fts.momentum().basicVector()) > 0);
00257 
00258 }
00259 
00260 
00261 //
00262 // calculate pseudorapidity from r and z
00263 //
00264 float DirectTrackerNavigation::calculateEta(float r, float z) const {
00265 
00266   if ( z > 0 ) return -log((tan(atan(r/z)/2.)));
00267   return log(-(tan(atan(r/z)/2.)));
00268 
00269 }