LDMX Software
EcalGeometry.h
Go to the documentation of this file.
1
15#ifndef DETDESCR_ECALGEOMETRY_H_
16#define DETDESCR_ECALGEOMETRY_H_
17
18// LDMX
19#include "DetDescr/EcalID.h"
21#include "Framework/Configure/Parameters.h"
22#include "Framework/Exception/Exception.h"
23
24// STL
25#include <map>
26
27// ROOT
28#include "TH2Poly.h"
29
30namespace ecal {
31class EcalGeometryProvider;
32}
33
34namespace ldmx {
35
102 public:
103 static constexpr const char* CONDITIONS_OBJECT_NAME{"EcalGeometry"};
104
111 virtual ~EcalGeometry() = default;
112
129 EcalID getID(double x, double y, double z, bool fallible = false) const;
130
150 EcalID getID(double x, double y, int layer_id, bool fallible = false) const;
151
175 EcalID getID(double x, double y, int layer_id, int module_id,
176 bool fallible = false) const;
177
187 std::tuple<double, double, double> getPosition(EcalID id) const;
188
192 std::pair<double, double> getPositionInModule(int cell_id) const;
193
200 double getZPosition(int layer) const {
201 return std::get<2>(layer_pos_xy_.at(layer));
202 }
203
209 int getNumLayers() const { return layer_pos_xy_.size(); }
210
223 int getNumModulesPerLayer() const { return 7; }
224
232 return cell_id_in_module_.GetNumberOfBins();
233 }
234
241 std::vector<EcalID> getNN(EcalID id) const {
242 auto list = NNMap_.at(EcalID(0, id.module(), id.cell()));
243 for (auto& flat : list)
244 flat = EcalID(id.layer(), flat.module(), flat.cell());
245 return list;
246 }
247
255 bool isNN(EcalID centroid, EcalID probe) const {
256 for (auto& id : getNN(centroid)) {
257 if (id == probe) return true;
258 }
259 return false;
260 }
261
268 std::vector<EcalID> getNNN(EcalID id) const {
269 auto list = NNNMap_.at(EcalID(0, id.module(), id.cell()));
270 for (auto& flat : list)
271 flat = EcalID(id.layer(), flat.module(), flat.cell());
272 return list;
273 }
274
283 bool isNNN(EcalID centroid, EcalID probe) const {
284 for (auto& id : getNNN(centroid)) {
285 if (id == probe) return true;
286 }
287 return false;
288 }
289
295 double getModuleMinR() const { return moduler_; }
296
302 double getModuleMaxR() const { return moduleR_; }
303
309 double getCellMinR() const { return cellr_; }
310
316 double getCellMaxR() const { return cellR_; }
317
327 TH2Poly* getCellPolyMap() const {
328 for (auto const& [cell_id, cell_center] : cell_pos_in_module_) {
329 cell_id_in_module_.Fill(cell_center.first, cell_center.second, cell_id);
330 }
331 return &cell_id_in_module_;
332 }
333
334 static EcalGeometry* debugMake(const framework::config::Parameters& p) {
335 return new EcalGeometry(p);
336 }
337
338 private:
344 EcalGeometry(const framework::config::Parameters& ps);
345 friend class ecal::EcalGeometryProvider;
346
350 void buildLayerMap();
351
369 void buildModuleMap();
370
403 void buildCellMap();
404
418 void buildCellModuleMap();
419
438 void buildNeighborMaps();
439
448 double distanceToEdge(EcalID id) const;
449
458 bool isEdgeCell(EcalID cellModuleID) const {
459 return (distanceToEdge(cellModuleID) < cellR_);
460 }
461
477 bool isInside(double normX, double normY) const;
478
479 private:
481 double gap_;
482
484 double cellr_{0};
485
487 double moduler_{0};
488
490 double cellR_{0};
491
493 double moduleR_{0};
494
501
506
511
520
534
542 double nCellRHeight_{0};
543
545 double ecalFrontZ_{0};
546
548 std::vector<double> layerZPositions_;
549
550 private:
553
558 std::map<int, std::tuple<double, double, double>> layer_pos_xy_;
559
566 std::map<int, std::pair<double, double>> module_pos_xy_;
567
574 std::map<int, std::pair<double, double>> cell_pos_in_module_;
575
585 std::map<EcalID, std::pair<double, double>> cell_pos_in_layer_;
586
596 std::map<EcalID, std::tuple<double, double, double>> cell_global_pos_;
597
603 std::map<EcalID, std::vector<EcalID>> NNMap_;
604
610 std::map<EcalID, std::vector<EcalID>> NNNMap_;
611
620 mutable TH2Poly cell_id_in_module_;
621};
622
623} // namespace ldmx
624
625#endif
Base class for conditions information like pedestals, gains, electronics maps, etc.
Class that defines an ECal detector ID with a cell number.
Base class for all conditions objects, very simple.
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:27
Translation between real-space positions and cell IDs within the ECal.
std::vector< EcalID > getNN(EcalID id) const
Get the Nearest Neighbors of the input ID.
std::map< EcalID, std::tuple< double, double, double > > cell_global_pos_
Position of cell centers relative to world geometry.
std::map< EcalID, std::vector< EcalID > > NNNMap_
Map of cell ID to neighbors of neighbor cells.
EcalID getID(double x, double y, double z, bool fallible=false) const
Get a cell's ID number from its position.
double getModuleMinR() const
Get the center-to-flat radius of the module hexagons.
std::tuple< double, double, double > getPosition(EcalID id) const
Get a cell's position from its ID number.
int getNumLayers() const
Get the number of layers in the Ecal Geometry.
double getZPosition(int layer) const
Get the z-coordinate given the layer id.
void buildCellModuleMap()
Constructs the positions of all the cells in a layer relative to the ecal center.
double layer_shift_y_
shift of layers in the y-direction [mm]
double distanceToEdge(EcalID id) const
Distance to module edge, and whether cell is on edge of module.
std::pair< double, double > getPositionInModule(int cell_id) const
Get a cell's position within a module.
std::map< int, std::pair< double, double > > module_pos_xy_
Postion of module centers relative to the center of the layer in world coordinates.
void buildCellMap()
Constructs the flat-bottomed hexagonal grid (cellID) of corner-down hexagonal cells.
double moduler_
Center-to-Flat Radius of module hexagon [mm].
double cellr_
Center-to-Flat Radius of cell hexagon [mm].
bool layer_shift_odd_
shift odd layers
int getNumModulesPerLayer() const
Get the number of modules in the Ecal flower.
std::vector< double > layerZPositions_
The layer Z postions are with respect to the front of the ECal [mm].
double nCellRHeight_
Number of cell center-to-corner radii (one side of the cell) from the bottom to the top of the module...
std::map< int, std::tuple< double, double, double > > layer_pos_xy_
Position of layer centers in world coordinates (uses layer ID as key)
int getNumCellsPerModule() const
Get the number of cells in each module of the Ecal Geometry.
double ecalFrontZ_
Front of ECal relative to world geometry [mm].
std::map< int, std::pair< double, double > > cell_pos_in_module_
Position of cell centers relative to center of module in p,q space.
double layer_shift_x_
shift of layers in the x-direction [mm]
std::vector< EcalID > getNNN(EcalID id) const
Get the Next-to-Nearest Neighbors of the input ID.
double cellR_
Center-to-Corner Radius of cell hexagon [mm].
double getModuleMaxR() const
Get the center-to-corner radius of the module hexagons.
void buildModuleMap()
Constructs the positions of the seven modules (moduleID) relative to the ecal center.
std::map< EcalID, std::pair< double, double > > cell_pos_in_layer_
Position of cell centers relative to center of layer in world coordinates.
bool isNNN(EcalID centroid, EcalID probe) const
Check if the probe id is one of the next-to-nearest neightbors of the centroid id.
int verbose_
verbosity
double getCellMinR() const
Get the center-to-flat radius of the cell hexagons.
bool isInside(double normX, double normY) const
Determines if point (x,y), already normed to max hexagon radius, lies within a hexagon.
bool layer_shift_odd_bilayer_
shift odd bi layers
bool cornersSideUp_
indicator of geometry orientation if true, flower shape's corners side (ie: side with two modules) is...
void buildLayerMap()
Constructs the positions of the layers in world coordinates.
std::map< EcalID, std::vector< EcalID > > NNMap_
Map of cell ID to neighboring cells.
double moduleR_
Center-to-Corner Radius of module hexagon [mm].
double gap_
Gap between module flat sides [mm].
virtual ~EcalGeometry()=default
Class destructor.
TH2Poly cell_id_in_module_
Honeycomb Binning from ROOT.
bool isNN(EcalID centroid, EcalID probe) const
Check if the probe id is one of the nearest neightbors of the centroid id.
double getCellMaxR() const
Get the center-to-corner radius of the cell hexagons.
bool isEdgeCell(EcalID cellModuleID) const
Check if input cell is on the edge of a module.
void buildNeighborMaps()
Construts NNMap and NNNMap.
TH2Poly * getCellPolyMap() const
Get a reference to the TH2Poly used for Cell IDs.
Extension of DetectorID providing access to ECal layers and cell numbers in a hex grid.
Definition EcalID.h:20