CMS 3D CMS Logo

TrackerTreeGenerator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: TrackerTreeGenerator
4 // Class: TrackerTreeGenerator
5 //
11 //
12 // Original Author: Johannes Hauk
13 // Created: Fri Jan 16 14:09:52 CET 2009
14 // Modified by: Gregor Mittag (DESY)
15 //
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
31 
35 
44 
55 
58 
59 #include "TTree.h"
60 //
61 // class decleration
62 //
63 
64 class TrackerTreeGenerator : public edm::one::EDAnalyzer<edm::one::SharedResources> {
65 public:
66  explicit TrackerTreeGenerator(const edm::ParameterSet&);
67  ~TrackerTreeGenerator() override = default;
68 
69 private:
70  void beginJob() override;
71  void analyze(const edm::Event&, const edm::EventSetup&) override;
72  void endJob() override;
73 
74  // ----------member data ---------------------------
78 
80  std::vector<TrackerTreeVariables> vTkTreeVar_;
82 };
83 
84 //
85 // constants, enums and typedefs
86 //
87 
88 //
89 // static data member definitions
90 //
91 
92 //
93 // constructors and destructor
94 //
96  : geomDetToken_(esConsumes()),
97  ptpToken_(esConsumes()),
98  topoToken_(esConsumes()),
99  createEntryForDoubleSidedModule_(config.getParameter<bool>("createEntryForDoubleSidedModule")),
100  config_(config) {
101  usesResource(TFileService::kSharedResource);
102 }
103 
104 //
105 // member functions
106 //
107 
108 // ------------ method called to for each event ------------
110  // now try to take directly the ideal geometry independent of used geometry in Global Tag
111  const GeometricDet* geometricDet = &iSetup.getData(geomDetToken_);
112  const PTrackerParameters& ptp = iSetup.getData(ptpToken_);
113  const TrackerTopology* tTopo = &iSetup.getData(topoToken_);
114 
115  TrackerGeomBuilderFromGeometricDet trackerBuilder;
116  const TrackerGeometry* tkGeom = trackerBuilder.build(geometricDet, ptp, tTopo);
117  AlignableTracker alignableTracker{tkGeom, tTopo};
118  const auto& ns = alignableTracker.trackerNameSpace();
119 
120  edm::LogInfo("TrackerTreeGenerator") << "@SUB=TrackerTreeGenerator::analyze"
121  << "There are " << tkGeom->detIds().size() << " dets and "
122  << tkGeom->detUnitIds().size() << " detUnits in the Geometry Record";
123 
125  edm::LogInfo("TrackerTreeGenerator") << "@SUB=TrackerTreeGenerator::analyze"
126  << "Create entry for each module AND one entry for virtual "
127  << "double-sided module in addition";
128  } else {
129  edm::LogInfo("TrackerTreeGenerator") << "@SUB=TrackerTreeGenerator::analyze"
130  << "Create one entry for each physical module, do NOT create additional "
131  << "entry for virtual double-sided module";
132  }
133 
134  for (const auto& detId : tkGeom->detIds()) {
135  const GeomDet& geomDet = *tkGeom->idToDet(detId);
136  const Surface& surface = geomDet.surface();
137 
138  TrackerTreeVariables tkTreeVar;
139  const auto rawId = detId.rawId();
140  tkTreeVar.rawId = rawId;
141  tkTreeVar.subdetId = detId.subdetId();
142 
143  switch (tkTreeVar.subdetId) {
145  tkTreeVar.layer = tTopo->pxbLayer(detId);
146  tkTreeVar.half = ns.tpb().halfBarrelNumber(rawId);
147  tkTreeVar.rod = tTopo->pxbLadder(detId); // ... so, ladder is not per halfBarrel-Layer, but per barrel-layer!
148  tkTreeVar.module = tTopo->pxbModule(detId);
149  break;
151  tkTreeVar.layer = tTopo->pxfDisk(detId);
152  tkTreeVar.side = tTopo->pxfSide(detId);
153  tkTreeVar.half = ns.tpe().halfCylinderNumber(rawId);
154  tkTreeVar.blade = tTopo->pxfBlade(detId);
155  tkTreeVar.panel = tTopo->pxfPanel(detId);
156  tkTreeVar.module = tTopo->pxfModule(detId);
157  break;
159  tkTreeVar.layer = tTopo->tibLayer(detId);
160  tkTreeVar.side = tTopo->tibStringInfo(detId)[0];
161  tkTreeVar.half = ns.tib().halfShellNumber(rawId);
162  tkTreeVar.rod = tTopo->tibStringInfo(detId)[2];
163  tkTreeVar.outerInner = tTopo->tibStringInfo(detId)[1];
164  tkTreeVar.module = tTopo->tibModule(detId);
165  tkTreeVar.isDoubleSide = tTopo->tibIsDoubleSide(detId);
166  tkTreeVar.isRPhi = tTopo->tibIsRPhi(detId);
167  tkTreeVar.isStereo = tTopo->tibIsStereo(detId);
168  break;
170  tkTreeVar.layer = tTopo->tidWheel(detId);
171  tkTreeVar.side = tTopo->tidSide(detId);
172  tkTreeVar.ring = tTopo->tidRing(detId);
173  tkTreeVar.outerInner = tTopo->tidModuleInfo(detId)[0];
174  tkTreeVar.module = tTopo->tidModuleInfo(detId)[1];
175  tkTreeVar.isDoubleSide = tTopo->tidIsDoubleSide(detId);
176  tkTreeVar.isRPhi = tTopo->tidIsRPhi(detId);
177  tkTreeVar.isStereo = tTopo->tidIsStereo(detId);
178  break;
180  tkTreeVar.layer = tTopo->tobLayer(detId);
181  tkTreeVar.side = tTopo->tobRodInfo(detId)[0];
182  tkTreeVar.rod = tTopo->tobRodInfo(detId)[1];
183  tkTreeVar.module = tTopo->tobModule(detId);
184  tkTreeVar.isDoubleSide = tTopo->tobIsDoubleSide(detId);
185  tkTreeVar.isRPhi = tTopo->tobIsRPhi(detId);
186  tkTreeVar.isStereo = tTopo->tobIsStereo(detId);
187  break;
189  tkTreeVar.layer = tTopo->tecWheel(detId);
190  tkTreeVar.side = tTopo->tecSide(detId);
191  tkTreeVar.ring = tTopo->tecRing(detId);
192  tkTreeVar.petal = tTopo->tecPetalInfo(detId)[1];
193  tkTreeVar.outerInner = tTopo->tecPetalInfo(detId)[0];
194  tkTreeVar.module = tTopo->tecModule(detId);
195  tkTreeVar.isDoubleSide = tTopo->tecIsDoubleSide(detId);
196  tkTreeVar.isRPhi = tTopo->tecIsRPhi(detId);
197  tkTreeVar.isStereo = tTopo->tecIsStereo(detId);
198  break;
199  }
200 
201  LocalPoint lPModule(0., 0., 0.), lUDirection(1., 0., 0.), lVDirection(0., 1., 0.), lWDirection(0., 0., 1.);
202  GlobalPoint gPModule = surface.toGlobal(lPModule), gUDirection = surface.toGlobal(lUDirection),
203  gVDirection = surface.toGlobal(lVDirection), gWDirection = surface.toGlobal(lWDirection);
204  double dR(999.), dPhi(999.), dZ(999.);
205  switch (tkTreeVar.subdetId) {
209  dR = gWDirection.perp() - gPModule.perp();
210  dPhi = deltaPhi(gUDirection.barePhi(), gPModule.barePhi());
211  dZ = gVDirection.z() - gPModule.z();
212  tkTreeVar.uDirection = dPhi > 0. ? 1 : -1;
213  tkTreeVar.vDirection = dZ > 0. ? 1 : -1;
214  tkTreeVar.wDirection = dR > 0. ? 1 : -1;
215  break;
217  dR = gUDirection.perp() - gPModule.perp();
218  dPhi = deltaPhi(gVDirection.barePhi(), gPModule.barePhi());
219  dZ = gWDirection.z() - gPModule.z();
220  tkTreeVar.uDirection = dR > 0. ? 1 : -1;
221  tkTreeVar.vDirection = dPhi > 0. ? 1 : -1;
222  tkTreeVar.wDirection = dZ > 0. ? 1 : -1;
223  break;
226  dR = gVDirection.perp() - gPModule.perp();
227  dPhi = deltaPhi(gUDirection.barePhi(), gPModule.barePhi());
228  dZ = gWDirection.z() - gPModule.z();
229  tkTreeVar.uDirection = dPhi > 0. ? 1 : -1;
230  tkTreeVar.vDirection = dR > 0. ? 1 : -1;
231  tkTreeVar.wDirection = dZ > 0. ? 1 : -1;
232  break;
233  }
234  tkTreeVar.posR = gPModule.perp();
235  tkTreeVar.posPhi = gPModule.barePhi(); // = gPModule.barePhi().degrees();
236  tkTreeVar.posEta = gPModule.eta();
237  tkTreeVar.posX = gPModule.x();
238  tkTreeVar.posY = gPModule.y();
239  tkTreeVar.posZ = gPModule.z();
240 
241  if (auto stripGeomDetUnit = dynamic_cast<const StripGeomDetUnit*>(&geomDet)) { //is it a single physical module?
242  switch (tkTreeVar.subdetId) {
247  auto& topol = dynamic_cast<const StripTopology&>(stripGeomDetUnit->specificTopology());
248  tkTreeVar.nStrips = topol.nstrips();
249  break;
250  }
251  }
252 
254  // do so only for individual modules and not also one entry for the combined doubleSided Module
255  if (tkTreeVar.isDoubleSide)
256  continue;
257  }
258  vTkTreeVar_.push_back(tkTreeVar);
259  }
260 }
261 
262 // ------------ method called once each job just before starting event loop ------------
264 
265 // ------------ method called once each job just after ending the event loop ------------
267  UInt_t rawId(999), subdetId(999), layer(999), side(999), half(999), rod(999), ring(999), petal(999), blade(999),
268  panel(999), outerInner(999), module(999), nStrips(999);
269  Bool_t isDoubleSide(false), isRPhi(false), isStereo(false);
270  Int_t uDirection(999), vDirection(999), wDirection(999);
271  Float_t posR(999.F), posPhi(999.F), posEta(999.F), posX(999.F), posY(999.F), posZ(999.F);
272 
273  edm::Service<TFileService> fileService;
274  TFileDirectory treeDir = fileService->mkdir("TrackerTree");
275  auto trackerTree{treeDir.make<TTree>("TrackerTree", "IDs of all modules (ideal geometry)")};
276  trackerTree->Branch("RawId", &rawId, "RawId/i");
277  trackerTree->Branch("SubdetId", &subdetId, "SubdetId/i");
278  trackerTree->Branch("Layer", &layer, "Layer/i"); // Barrel: Layer, Forward: Disk
279  trackerTree->Branch("Side", &side, "Side/i"); // Rod/Ring in +z or -z
280  trackerTree->Branch("Half", &half, "Half/i"); // PXB: HalfBarrel, PXF: HalfCylinder, TIB: HalfShell
281  trackerTree->Branch("Rod", &rod, "Rod/i"); // Barrel (Ladder or String or Rod)
282  trackerTree->Branch("Ring", &ring, "Ring/i"); // Forward
283  trackerTree->Branch("Petal", &petal, "Petal/i"); // TEC
284  trackerTree->Branch("Blade", &blade, "Blade/i"); // PXF
285  trackerTree->Branch("Panel", &panel, "Panel/i"); // PXF
286  trackerTree->Branch("OuterInner", &outerInner, "OuterInner/i"); // front/back String,Ring,Petal
287  trackerTree->Branch("Module", &module, "Module/i"); // Module ID
288  trackerTree->Branch("NStrips", &nStrips, "NStrips/i");
289  trackerTree->Branch("IsDoubleSide", &isDoubleSide, "IsDoubleSide/O");
290  trackerTree->Branch("IsRPhi", &isRPhi, "IsRPhi/O");
291  trackerTree->Branch("IsStereo", &isStereo, "IsStereo/O");
292  trackerTree->Branch("UDirection", &uDirection, "UDirection/I");
293  trackerTree->Branch("VDirection", &vDirection, "VDirection/I");
294  trackerTree->Branch("WDirection", &wDirection, "WDirection/I");
295  trackerTree->Branch("PosR", &posR, "PosR/F");
296  trackerTree->Branch("PosPhi", &posPhi, "PosPhi/F");
297  trackerTree->Branch("PosEta", &posEta, "PosEta/F");
298  trackerTree->Branch("PosX", &posX, "PosX/F");
299  trackerTree->Branch("PosY", &posY, "PosY/F");
300  trackerTree->Branch("PosZ", &posZ, "PosZ/F");
301 
302  for (const auto& iTree : vTkTreeVar_) {
303  rawId = iTree.rawId;
304  subdetId = iTree.subdetId;
305  layer = iTree.layer;
306  side = iTree.side;
307  half = iTree.half;
308  rod = iTree.rod;
309  ring = iTree.ring;
310  petal = iTree.petal;
311  blade = iTree.blade;
312  panel = iTree.panel;
313  outerInner = iTree.outerInner;
314  module = iTree.module;
315  nStrips = iTree.nStrips;
316  isDoubleSide = iTree.isDoubleSide;
317  isRPhi = iTree.isRPhi;
318  isStereo = iTree.isStereo;
319  uDirection = iTree.uDirection;
320  vDirection = iTree.vDirection;
321  wDirection = iTree.wDirection;
322  posR = iTree.posR;
323  posPhi = iTree.posPhi;
324  posEta = iTree.posEta;
325  posX = iTree.posX;
326  posY = iTree.posY;
327  posZ = iTree.posZ;
328 
329  trackerTree->Fill();
330  }
331  edm::LogInfo("TrackerTreeGenerator") << "@SUB=TrackerTreeGenerator::endJob"
332  << "TrackerTree contains " << vTkTreeVar_.size() << " entries overall";
333 }
334 
335 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
static constexpr auto TEC
bool tibIsDoubleSide(const DetId &id) const
bool tecIsDoubleSide(const DetId &id) const
bool tidIsDoubleSide(const DetId &id) const
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
unsigned int tobLayer(const DetId &id) const
unsigned int pxbLayer(const DetId &id) const
T perp() const
Definition: PV3DBase.h:69
const edm::ESGetToken< GeometricDet, IdealGeometryRecord > geomDetToken_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
TrackerTreeGenerator(const edm::ParameterSet &)
~TrackerTreeGenerator() override=default
std::vector< unsigned int > tidModuleInfo(const DetId &id) const
unsigned int pxfBlade(const DetId &id) const
T z() const
Definition: PV3DBase.h:61
const DetIdContainer & detIds() const override
Returm a vector of all GeomDet DetIds (including those of GeomDetUnits)
unsigned int tibModule(const DetId &id) const
unsigned int tidSide(const DetId &id) const
unsigned int pxfModule(const DetId &id) const
std::vector< unsigned int > tecPetalInfo(const DetId &id) const
unsigned int tidWheel(const DetId &id) const
unsigned int tecWheel(const DetId &id) const
T eta() const
Definition: PV3DBase.h:73
TrackerGeometry * build(const GeometricDet *gd, const PTrackerParameters &ptp, const TrackerTopology *tTopo)
Definition: config.py:1
unsigned int pxbLadder(const DetId &id) const
edm::ParameterSet config_
void analyze(const edm::Event &, const edm::EventSetup &) override
unsigned int tecRing(const DetId &id) const
ring id
bool tibIsStereo(const DetId &id) const
bool tobIsDoubleSide(const DetId &id) const
nStrips
1.2 is to make the matching window safely the two nearest strips 0.35 is the size of an ME0 chamber i...
T barePhi() const
Definition: PV3DBase.h:65
bool tobIsStereo(const DetId &id) const
bool tibIsRPhi(const DetId &id) const
T * make(const Args &...args) const
make new ROOT object
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
unsigned int tecModule(const DetId &id) const
int iEvent
Definition: GenABIO.cc:224
unsigned int tecSide(const DetId &id) const
const edm::ESGetToken< PTrackerParameters, PTrackerParametersRcd > ptpToken_
unsigned int pxfDisk(const DetId &id) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static constexpr auto TOB
const TrackerGeomDet * idToDet(DetId) const override
bool tobIsRPhi(const DetId &id) const
bool tidIsRPhi(const DetId &id) const
bool tecIsRPhi(const DetId &id) const
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
unsigned int pxfPanel(const DetId &id) const
Log< level::Info, false > LogInfo
unsigned int pxfSide(const DetId &id) const
std::vector< unsigned int > tibStringInfo(const DetId &id) const
static constexpr auto TIB
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
bool tecIsStereo(const DetId &id) const
bool tidIsStereo(const DetId &id) const
std::vector< unsigned int > tobRodInfo(const DetId &id) const
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
const bool createEntryForDoubleSidedModule_
std::vector< TrackerTreeVariables > vTkTreeVar_
unsigned int tidRing(const DetId &id) const
unsigned int tibLayer(const DetId &id) const
unsigned int tobModule(const DetId &id) const
unsigned int pxbModule(const DetId &id) const
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
static constexpr auto TID
const DetIdContainer & detUnitIds() const override
Returm a vector of all GeomDetUnit DetIds.