File indexing completed on 2024-04-06 12:26:55
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "RecoMuon/GlobalTrackingTools/interface/DirectTrackerNavigation.h"
0010
0011
0012
0013
0014
0015 #include <algorithm>
0016
0017
0018
0019
0020
0021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0022 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
0023 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
0024 #include "DataFormats/GeometrySurface/interface/BoundCylinder.h"
0025 #include "DataFormats/GeometrySurface/interface/BoundDisk.h"
0026
0027 using namespace std;
0028
0029
0030
0031
0032 DirectTrackerNavigation::DirectTrackerNavigation(const edm::ESHandle<GeometricSearchTracker>& tkLayout, bool outOnly)
0033 : theGeometricSearchTracker(tkLayout), theOutLayerOnlyFlag(outOnly), theEpsilon(-0.01) {}
0034
0035
0036
0037
0038 vector<const DetLayer*> DirectTrackerNavigation::compatibleLayers(const FreeTrajectoryState& fts,
0039 PropagationDirection dir) const {
0040 bool inOut = outward(fts);
0041 double eta0 = fts.position().eta();
0042
0043 vector<const DetLayer*> output;
0044
0045
0046
0047 if (inOut) {
0048 if (!theOutLayerOnlyFlag) {
0049 inOutPx(fts, output);
0050
0051 if (eta0 > 1.55)
0052 inOutFPx(fts, output);
0053 else if (eta0 < -1.55)
0054 inOutBPx(fts, output);
0055
0056 if (fabs(eta0) < 1.67)
0057 inOutTIB(fts, output);
0058
0059 if (eta0 > 1.17)
0060 inOutFTID(fts, output);
0061 else if (eta0 < -1.17)
0062 inOutBTID(fts, output);
0063 }
0064
0065 if (fabs(eta0) < 1.35)
0066 inOutTOB(fts, output);
0067
0068 if (eta0 > 0.97)
0069 inOutFTEC(fts, output);
0070 else if (eta0 < -0.97)
0071 inOutBTEC(fts, output);
0072
0073 } else {
0074 LogTrace("Muon|RecoMuon|DirectionTrackerNavigation") << "No implementation for inward state at this moment. ";
0075 }
0076
0077 if (dir == oppositeToMomentum)
0078 std::reverse(output.begin(), output.end());
0079
0080 return output;
0081 }
0082
0083
0084
0085
0086 void DirectTrackerNavigation::inOutPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0087 for (const auto i : theGeometricSearchTracker->pixelBarrelLayers()) {
0088 if (checkCompatible(fts, i))
0089 output.push_back(i);
0090 }
0091 }
0092
0093
0094
0095
0096 void DirectTrackerNavigation::inOutTIB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0097 for (const auto i : theGeometricSearchTracker->tibLayers()) {
0098 if (checkCompatible(fts, i))
0099 output.push_back(i);
0100 }
0101 }
0102
0103
0104
0105
0106 void DirectTrackerNavigation::inOutTOB(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0107 for (const auto i : theGeometricSearchTracker->tobLayers()) {
0108 if (checkCompatible(fts, i))
0109 output.push_back(i);
0110 }
0111 }
0112
0113
0114
0115
0116 void DirectTrackerNavigation::inOutFPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0117 for (const auto i : theGeometricSearchTracker->posPixelForwardLayers()) {
0118 if (checkCompatible(fts, i))
0119 output.push_back(i);
0120 }
0121 }
0122
0123
0124
0125
0126 void DirectTrackerNavigation::inOutFTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0127 for (const auto i : theGeometricSearchTracker->posTidLayers()) {
0128 if (checkCompatible(fts, i))
0129 output.push_back(i);
0130 }
0131 }
0132
0133
0134
0135
0136 void DirectTrackerNavigation::inOutFTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0137 for (const auto i : theGeometricSearchTracker->posTecLayers()) {
0138 if (checkCompatible(fts, i))
0139 output.push_back(i);
0140 }
0141 }
0142
0143
0144
0145
0146 void DirectTrackerNavigation::inOutBPx(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0147 for (const auto i : theGeometricSearchTracker->negPixelForwardLayers()) {
0148 if (checkCompatible(fts, i))
0149 output.push_back(i);
0150 }
0151 }
0152
0153
0154
0155
0156 void DirectTrackerNavigation::inOutBTID(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0157 for (const auto i : theGeometricSearchTracker->negTidLayers()) {
0158 if (checkCompatible(fts, i))
0159 output.push_back(i);
0160 }
0161 }
0162
0163
0164
0165
0166 void DirectTrackerNavigation::inOutBTEC(const FreeTrajectoryState& fts, vector<const DetLayer*>& output) const {
0167 for (const auto i : theGeometricSearchTracker->negTecLayers()) {
0168 if (checkCompatible(fts, i))
0169 output.push_back(i);
0170 }
0171 }
0172
0173
0174
0175
0176 bool DirectTrackerNavigation::checkCompatible(const FreeTrajectoryState& fts, const BarrelDetLayer* dl) const {
0177 float eta0 = fts.position().eta();
0178
0179 const BoundCylinder& bc = dl->specificSurface();
0180 float radius = bc.radius();
0181 float length = bc.bounds().length() / 2.;
0182
0183 float eta = calculateEta(radius, length);
0184
0185 return (fabs(eta0) <= (fabs(eta) + theEpsilon));
0186 }
0187
0188
0189
0190
0191 bool DirectTrackerNavigation::checkCompatible(const FreeTrajectoryState& fts, const ForwardDetLayer* dl) const {
0192 float eta0 = fts.position().eta();
0193
0194 const BoundDisk& bd = dl->specificSurface();
0195
0196 float outRadius = bd.outerRadius();
0197 float inRadius = bd.innerRadius();
0198 float z = bd.position().z();
0199
0200 float etaOut = calculateEta(outRadius, z);
0201 float etaIn = calculateEta(inRadius, z);
0202
0203 if (eta0 > 0)
0204 return (eta0 > (etaOut - theEpsilon) && eta0 < (etaIn + theEpsilon));
0205 else
0206 return (eta0 < (etaOut + theEpsilon) && eta0 > (etaIn - theEpsilon));
0207 }
0208
0209
0210
0211
0212 bool DirectTrackerNavigation::outward(const FreeTrajectoryState& fts) const {
0213 return (fts.position().basicVector().dot(fts.momentum().basicVector()) > 0);
0214 }
0215
0216
0217
0218
0219 float DirectTrackerNavigation::calculateEta(float r, float z) const {
0220 if (z > 0)
0221 return -log((tan(atan(r / z) / 2.)));
0222 return log(-(tan(atan(r / z) / 2.)));
0223 }