CMS 3D CMS Logo

GeometricSearchTracker.cc
Go to the documentation of this file.
3 
5 
6 using namespace std;
7 
8 GeometricSearchTracker::GeometricSearchTracker(const vector<BarrelDetLayer const*>& pxlBar,
9  const vector<BarrelDetLayer const*>& tib,
10  const vector<BarrelDetLayer const*>& tob,
11  const vector<ForwardDetLayer const*>& negPxlFwd,
12  const vector<ForwardDetLayer const*>& negTid,
13  const vector<ForwardDetLayer const*>& negTec,
14  const vector<ForwardDetLayer const*>& posPxlFwd,
15  const vector<ForwardDetLayer const*>& posTid,
16  const vector<ForwardDetLayer const*>& posTec,
17  const TrackerTopology* tTopo)
18  : thePixelBarrelLayers(pxlBar.begin(), pxlBar.end()),
19  theTibLayers(tib.begin(), tib.end()),
20  theTobLayers(tob.begin(), tob.end()),
21  theNegPixelForwardLayers(negPxlFwd.begin(), negPxlFwd.end()),
22  theNegTidLayers(negTid.begin(), negTid.end()),
23  theNegTecLayers(negTec.begin(), negTec.end()),
24  thePosPixelForwardLayers(posPxlFwd.begin(), posPxlFwd.end()),
25  thePosTidLayers(posTid.begin(), posTid.end()),
26  thePosTecLayers(posTec.begin(), posTec.end()),
27  theTrkTopo(tTopo) {
29  theBarrelLayers.insert(theBarrelLayers.end(), theTibLayers.begin(), theTibLayers.end());
30  theBarrelLayers.insert(theBarrelLayers.end(), theTobLayers.begin(), theTobLayers.end());
31 
35 
39 
42  theAllLayers.assign(theBarrelLayers.begin(), theBarrelLayers.end());
43  theAllLayers.insert(theAllLayers.end(), theForwardLayers.begin(), theForwardLayers.end());
44 
45  // number the layers
46  int sq = 0;
47  for (auto l : theAllLayers)
48  const_cast<DetLayer&>(*l).setSeqNum(sq++);
49 
50  edm::LogInfo("TkDetLayers") << "------ GeometricSearchTracker constructed with: ------"
51  << "\n"
52  << "n pxlBarLayers: " << this->pixelBarrelLayers().size() << "\n"
53  << "n tibLayers: " << this->tibLayers().size() << "\n"
54  << "n tobLayers: " << this->tobLayers().size() << "\n"
55  << "n negPxlFwdLayers: " << this->negPixelForwardLayers().size() << "\n"
56  << "n posPxlFwdLayers: " << this->posPixelForwardLayers().size() << "\n"
57  << "n negTidLayers: " << this->negTidLayers().size() << "\n"
58  << "n posTidLayers: " << this->posTidLayers().size() << "\n"
59  << "n negTecLayers: " << this->negTecLayers().size() << "\n"
60  << "n posTecLayers: " << this->posTecLayers().size() << "\n"
61 
62  << "n barreLayers: " << this->barrelLayers().size() << "\n"
63  << "n negforwardLayers: " << this->negForwardLayers().size() << "\n"
64  << "n posForwardLayers: " << this->posForwardLayers().size()
65  << "\nn Total : " << theAllLayers.size() << " " << sq << std::endl;
66 
67  for (auto l : theAllLayers)
68  edm::LogInfo("TkDetLayers") << (*l).seqNum() << ": " << (*l).subDetector() << ", ";
69  edm::LogInfo("TkDetLayers") << std::endl;
70 }
71 
73  for (auto l : theAllLayers)
74  delete const_cast<DetLayer*>(l);
75 }
76 
79 
81  //If it's a tracker Det
82  if (id.det() == DetId::Detector::Tracker) {
83  switch (id.subdetId()) {
85  return theTibLayers[theTrkTopo->tibLayer(id) - 1];
86  break;
87 
89  return theTobLayers[theTrkTopo->tobLayer(id) - 1];
90  break;
91 
93  if (theTrkTopo->tidSide(id) == 1) {
94  return theNegTidLayers[theTrkTopo->tidWheel(id) - 1];
95  } else if (theTrkTopo->tidSide(id) == 2) {
96  return thePosTidLayers[theTrkTopo->tidWheel(id) - 1];
97  }
98  break;
99 
101  if (theTrkTopo->tecSide(id) == 1) {
102  return theNegTecLayers[theTrkTopo->tecWheel(id) - 1];
103  } else if (theTrkTopo->tecSide(id) == 2) {
104  return thePosTecLayers[theTrkTopo->tecWheel(id) - 1];
105  }
106  break;
107 
109  return thePixelBarrelLayers[theTrkTopo->pxbLayer(id) - 1];
110  break;
111 
113  if (theTrkTopo->pxfSide(id) == 1) {
115  } else if (theTrkTopo->pxfSide(id) == 2) {
117  }
118  break;
119 
120  default:
121  edm::LogError("TkDetLayers") << "ERROR:layer not found!";
122  // throw(something);
123  }
124  return nullptr; //just to avoid compile warnings
125  } else if (id.det() == DetId::Forward && id.subdetId() == FastTime) {
126  //If it's MTD
127  return mtdDetLayerGeometry->idToLayer(id);
128  }
129  return nullptr; //just to avoid compile warnings
130 }
131 
133 
134 void GeometricSearchTracker::addMTDLayers(const std::vector<BarrelDetLayer const*>& btl,
135  const std::vector<ForwardDetLayer const*>& negEtl,
136  const std::vector<ForwardDetLayer const*>& posEtl) {
137  //Barrel
138  theBTLLayers.assign(btl.begin(), btl.end());
139  theBarrelLayers.insert(theBarrelLayers.end(), theBTLLayers.begin(), theBTLLayers.end());
140  //Endcap
141  theNegETLLayers.assign(negEtl.begin(), negEtl.end());
142  thePosETLLayers.assign(posEtl.begin(), posEtl.end());
143  theETLLayers.assign(negEtl.begin(), negEtl.end());
144  theETLLayers.insert(theETLLayers.end(), posEtl.begin(), posEtl.end());
147  //Reordering of tracker + MTD layers
148  theForwardLayers.clear();
149  theAllLayers.clear();
152  theAllLayers.assign(theBarrelLayers.begin(), theBarrelLayers.end());
153  theAllLayers.insert(theAllLayers.end(), theForwardLayers.begin(), theForwardLayers.end());
154 
155  // number the layers
156  int sq = 0;
157  for (auto l : theAllLayers)
158  const_cast<DetLayer&>(*l).setSeqNum(sq++);
159 
160  edm::LogInfo("MTDDetLayers") << "------ GeometricSearchTracker+MTD constructed with: ------"
161  << "\n"
162  << "n pxlBarLayers: " << this->pixelBarrelLayers().size() << "\n"
163  << "n tibLayers: " << this->tibLayers().size() << "\n"
164  << "n tobLayers: " << this->tobLayers().size() << "\n"
165  << "n negPxlFwdLayers: " << this->negPixelForwardLayers().size() << "\n"
166  << "n posPxlFwdLayers: " << this->posPixelForwardLayers().size() << "\n"
167  << "n negTidLayers: " << this->negTidLayers().size() << "\n"
168  << "n posTidLayers: " << this->posTidLayers().size() << "\n"
169  << "n negTecLayers: " << this->negTecLayers().size() << "\n"
170  << "n posTecLayers: " << this->posTecLayers().size() << "\n"
171  << "n barreLayers: " << this->barrelLayers().size() << "\n"
172  << "n negforwardLayers: " << this->negForwardLayers().size() << "\n"
173  << "n posForwardLayers: " << this->posForwardLayers().size() << "\n"
174  << "n MTDbarrelLayers: " << this->theBTLLayers.size() << "\n"
175  << "n MTDnegLayers: " << this->theNegETLLayers.size() << "\n"
176  << "n MTDposLayers: " << this->thePosETLLayers.size() << "\n"
177  << "\nn Total : " << theAllLayers.size() << std::endl;
178 }
PixelSubdetector.h
MessageLogger.h
PixelSubdetector::PixelEndcap
Definition: PixelSubdetector.h:11
GeometricSearchTracker::tobLayers
std::vector< BarrelDetLayer const * > const & tobLayers() const
Definition: GeometricSearchTracker.h:46
PixelSubdetector::PixelBarrel
Definition: PixelSubdetector.h:11
GeometricSearchTracker::tibLayers
std::vector< BarrelDetLayer const * > const & tibLayers() const
Definition: GeometricSearchTracker.h:45
GeometricSearchTracker::mtdDetLayerGeometry
MTDDetLayerGeometry * mtdDetLayerGeometry
Definition: GeometricSearchTracker.h:60
DetLayer
Definition: DetLayer.h:21
GeometricSearchTracker::addDetLayerGeometry
void addDetLayerGeometry()
Definition: GeometricSearchTracker.cc:132
TrackerTopology::pxfSide
unsigned int pxfSide(const DetId &id) const
Definition: TrackerTopology.h:192
GeometricSearchTracker::theTrkTopo
const TrackerTopology * theTrkTopo
Definition: GeometricSearchTracker.h:89
TrackerTopology
Definition: TrackerTopology.h:16
GeometricSearchTracker::GeometricSearchTracker
GeometricSearchTracker(const std::vector< BarrelDetLayer const * > &pxlBar, const std::vector< BarrelDetLayer const * > &tib, const std::vector< BarrelDetLayer const * > &tob, const std::vector< ForwardDetLayer const * > &negPxlFwd, const std::vector< ForwardDetLayer const * > &negTid, const std::vector< ForwardDetLayer const * > &negTec, const std::vector< ForwardDetLayer const * > &posPxlFwd, const std::vector< ForwardDetLayer const * > &posTid, const std::vector< ForwardDetLayer const * > &posTec, const TrackerTopology *tTopo) __attribute__((cold))
Definition: GeometricSearchTracker.cc:8
GeometricSearchTracker::thePosPixelForwardLayers
std::vector< ForwardDetLayer const * > thePosPixelForwardLayers
Definition: GeometricSearchTracker.h:79
GeometricSearchTracker::theTobLayers
std::vector< BarrelDetLayer const * > theTobLayers
Definition: GeometricSearchTracker.h:74
GeometricSearchTracker::thePixelBarrelLayers
std::vector< BarrelDetLayer const * > thePixelBarrelLayers
Definition: GeometricSearchTracker.h:72
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
GeometricSearchTracker::thePosTecLayers
std::vector< ForwardDetLayer const * > thePosTecLayers
Definition: GeometricSearchTracker.h:81
TrackerTopology::tidWheel
unsigned int tidWheel(const DetId &id) const
Definition: TrackerTopology.h:201
TrackerTopology::pxbLayer
unsigned int pxbLayer(const DetId &id) const
Definition: TrackerTopology.h:144
align::Tracker
Definition: StructureType.h:70
DetId
Definition: DetId.h:17
GeometricSearchTracker.h
GeometricSearchTracker::theNegTecLayers
std::vector< ForwardDetLayer const * > theNegTecLayers
Definition: GeometricSearchTracker.h:78
GeometricSearchTracker::thePosForwardLayers
std::vector< ForwardDetLayer const * > thePosForwardLayers
Definition: GeometricSearchTracker.h:70
TrackerTopology.h
GeometricSearchTracker::posTecLayers
std::vector< ForwardDetLayer const * > const & posTecLayers() const
Definition: GeometricSearchTracker.h:54
GeometricSearchTracker::thePosETLLayers
std::vector< ForwardDetLayer const * > thePosETLLayers
Definition: GeometricSearchTracker.h:87
GeometricSearchTracker::theNegTidLayers
std::vector< ForwardDetLayer const * > theNegTidLayers
Definition: GeometricSearchTracker.h:77
mps_fire.end
end
Definition: mps_fire.py:242
GeometricSearchTracker::~GeometricSearchTracker
~GeometricSearchTracker() override __attribute__((cold))
Definition: GeometricSearchTracker.cc:72
GeometricSearchTracker::posTidLayers
std::vector< ForwardDetLayer const * > const & posTidLayers() const
Definition: GeometricSearchTracker.h:53
GeometricSearchTracker::barrelLayers
std::vector< BarrelDetLayer const * > const & barrelLayers() const
Definition: GeometricSearchTracker.h:38
StripSubdetector::TIB
static constexpr auto TIB
Definition: StripSubdetector.h:16
GeometricSearchTracker::theNegPixelForwardLayers
std::vector< ForwardDetLayer const * > theNegPixelForwardLayers
Definition: GeometricSearchTracker.h:76
GeometricSearchTracker::thePosTidLayers
std::vector< ForwardDetLayer const * > thePosTidLayers
Definition: GeometricSearchTracker.h:80
GeometricSearchTracker::negTecLayers
std::vector< ForwardDetLayer const * > const & negTecLayers() const
Definition: GeometricSearchTracker.h:50
GeometricSearchTracker::theNegETLLayers
std::vector< ForwardDetLayer const * > theNegETLLayers
Definition: GeometricSearchTracker.h:86
GeometricSearchTracker::posPixelForwardLayers
std::vector< ForwardDetLayer const * > const & posPixelForwardLayers() const
Definition: GeometricSearchTracker.h:52
GeometricSearchTracker::theBTLLayers
std::vector< BarrelDetLayer const * > theBTLLayers
Definition: GeometricSearchTracker.h:84
GeometricSearchTracker::theETLLayers
std::vector< ForwardDetLayer const * > theETLLayers
Definition: GeometricSearchTracker.h:85
GeometricSearchTracker::theAllLayers
std::vector< DetLayer const * > theAllLayers
Definition: GeometricSearchTracker.h:66
FastTime
Definition: ForwardSubdetector.h:6
TrackerTopology::tidSide
unsigned int tidSide(const DetId &id) const
Definition: TrackerTopology.h:190
MTDDetLayerGeometry::idToLayer
const DetLayer * idToLayer(const DetId &detId) const override
return the DetLayer which correspond to a certain DetId
Definition: MTDDetLayerGeometry.cc:97
TrackerTopology::pxfDisk
unsigned int pxfDisk(const DetId &id) const
Definition: TrackerTopology.h:446
GeometricSearchTracker::theBarrelLayers
std::vector< BarrelDetLayer const * > theBarrelLayers
Definition: GeometricSearchTracker.h:67
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
MTDDetLayerGeometry
Definition: MTDDetLayerGeometry.h:22
cmsLHEtoEOSManager.l
l
Definition: cmsLHEtoEOSManager.py:204
GeometricSearchTracker::idToLayer
const DetLayer * idToLayer(const DetId &detId) const override
Give the DetId of a module, returns the pointer to the corresponding DetLayer.
Definition: GeometricSearchTracker.cc:80
TrackerTopology::tobLayer
unsigned int tobLayer(const DetId &id) const
Definition: TrackerTopology.h:147
GeometricSearchTracker::negForwardLayers
std::vector< ForwardDetLayer const * > const & negForwardLayers() const
Definition: GeometricSearchTracker.h:41
std
Definition: JetResolutionObject.h:76
StripSubdetector::TEC
static constexpr auto TEC
Definition: StripSubdetector.h:19
GeometricSearchTracker::theTibLayers
std::vector< BarrelDetLayer const * > theTibLayers
Definition: GeometricSearchTracker.h:73
GeometricSearchTracker::negPixelForwardLayers
std::vector< ForwardDetLayer const * > const & negPixelForwardLayers() const
Definition: GeometricSearchTracker.h:48
GeometricSearchTracker::addMTDLayers
void addMTDLayers(const std::vector< BarrelDetLayer const * > &btl, const std::vector< ForwardDetLayer const * > &negEtl, const std::vector< ForwardDetLayer const * > &posEtl)
Definition: GeometricSearchTracker.cc:134
GeometricSearchTracker::posForwardLayers
std::vector< ForwardDetLayer const * > const & posForwardLayers() const
Definition: GeometricSearchTracker.h:42
StripSubdetector::TOB
static constexpr auto TOB
Definition: StripSubdetector.h:18
GeometricSearchTracker::theForwardLayers
std::vector< ForwardDetLayer const * > theForwardLayers
Definition: GeometricSearchTracker.h:68
GeometricSearchTracker::theNegForwardLayers
std::vector< ForwardDetLayer const * > theNegForwardLayers
Definition: GeometricSearchTracker.h:69
DetId::Forward
Definition: DetId.h:30
TrackerTopology::tecWheel
unsigned int tecWheel(const DetId &id) const
Definition: TrackerTopology.h:198
StripSubdetector.h
StripSubdetector::TID
static constexpr auto TID
Definition: StripSubdetector.h:17
GeometricSearchTracker::negTidLayers
std::vector< ForwardDetLayer const * > const & negTidLayers() const
Definition: GeometricSearchTracker.h:49
GeometricSearchTracker::pixelBarrelLayers
std::vector< BarrelDetLayer const * > const & pixelBarrelLayers() const
Definition: GeometricSearchTracker.h:44
TrackerTopology::tibLayer
unsigned int tibLayer(const DetId &id) const
Definition: TrackerTopology.h:150
TrackerTopology::tecSide
unsigned int tecSide(const DetId &id) const
Definition: TrackerTopology.h:184