5static double distance(
const std::pair<double, double>& p1,
6 const std::pair<double, double>& p2) {
7 return sqrt((p1.first - p2.first) * (p1.first - p2.first) +
8 (p1.second - p2.second) * (p1.second - p2.second));
11[[maybe_unused]]
static double distance(
12 const std::tuple<double, double, double>& p1,
13 const std::tuple<double, double, double>& p2) {
14 return sqrt((std::get<0>(p1) - std::get<0>(p2)) *
15 (std::get<0>(p1) - std::get<0>(p2)) +
16 (std::get<1>(p1) - std::get<1>(p2)) *
17 (std::get<1>(p1) - std::get<1>(p2)) +
18 (std::get<2>(p1) - std::get<2>(p2)) *
19 (std::get<2>(p1) - std::get<2>(p2)));
25static void rotate(
double& p,
double& q) {
34static void unrotate(
double& p,
double& q) {
55 EXCEPTION_RAISE(
"BadConf",
56 "Cannot shift both odd sensitive layers and odd bilayers");
63 ldmx_log(debug) <<
"Building module map with gap " << std::setprecision(2)
65 <<
", min/max radii of cell " <<
cell_r_min_ <<
" / "
75 ldmx_log(trace) <<
"Geometry fully constructed";
91 EXCEPTION_RAISE(
"BadConf",
"z = " + std::to_string(z) +
92 " mm is not within any"
93 " of the configured layers.");
96 return getID(x, y, layer_id, fallible);
100 bool fallible)
const {
113 double probe_x{p - module_xy.first}, probe_y{q - module_xy.second};
128 "Coordinates relative to layer (p,q) = (%.2f, %.2f) mm "
129 "derived from world coordinates (%.2f, %.2f) mm with layer = %d "
130 "are not inside any module.",
131 p, q, x, y, layer_id)
136 return getID(x, y, layer_id, module_id);
140 bool fallible)
const {
161 "Relative cell coordinates (%.2f, %.2f) mm "
162 "derived from world coordinates (%.2f, %.2f) mm with layer = %d "
163 "and module = %d are outside module hexagon",
164 p, q, x, y, layer_id, module_id)
169 return EcalID(layer_id, module_id, cell_id);
190 ldmx_log(debug) <<
"Shifting odd layers by (x,y) = (" <<
layer_shift_x_
193 ldmx_log(debug) <<
"Shifting odd bilayers by (x,y) = (" <<
layer_shift_x_
208 ldmx_log(trace) <<
" Layer " << i_layer <<
" has center at (" << x
209 <<
", " << y <<
", " << z <<
") mm";
216 static const double c_pi = 3.14159265358979323846;
227 for (
unsigned id = 1;
id < 7;
id++) {
237 ldmx_log(trace) <<
" Module " <<
id <<
" is centered at (x,y) = " <<
"("
238 << x <<
", " << y <<
") mm";
265 double grid_min_p = 0., grid_min_q = 0.;
266 int num_p_cells = 0, num_q_cells = 0;
280 if (num_q_cells % 2 == 0)
290 grid_map.Honeycomb(grid_min_p, grid_min_q,
cell_r_max_, num_p_cells,
293 ldmx_log(trace) << std::setprecision(2)
294 <<
"Building buildCellMap with cell rmin: " <<
cell_r_min_
295 <<
" cell rmax: " <<
cell_r_max_ <<
" (gridMinP,gridMinQ) = ("
296 << grid_min_p <<
"," << grid_min_q <<
")"
297 <<
" (numPCells,numQCells) = (" << num_p_cells <<
","
298 << num_q_cells <<
")";
301 TListIter next(grid_map.GetBins());
302 TH2PolyBin* poly_bin = 0;
306 while ((poly_bin = (TH2PolyBin*)next())) {
310 poly = (TGraph*)poly_bin->GetPolygon();
314 int num_vertices_inside = 0;
315 double vertex_p[6], vertex_q[6];
317 ldmx_log(trace) <<
" Cell vertices";
318 for (
unsigned i = 0; i < 6; i++) {
319 poly->GetPoint(i, vertex_p[i], vertex_q[i]);
320 ldmx_log(trace) <<
" vtx # " << i;
321 ldmx_log(trace) <<
" vtx p,q " << vertex_p[i] <<
" " << vertex_q[i];
325 if (isinside[i]) num_vertices_inside++;
328 if (num_vertices_inside > 1) {
331 double actual_p[8], actual_q[8];
333 if (num_vertices_inside < 6) {
337 ldmx_log(trace) <<
" Polygon " << cell_id
338 <<
" has vertices poking out of module hexagon.";
341 for (
int i = 0; i < 6; i++) {
342 int up = i == 5 ? 0 : i + 1;
343 int dn = i == 0 ? 5 : i - 1;
344 if (isinside[i] and (not isinside[up] or not isinside[dn])) {
351 double edge_origin_p, edge_origin_q;
352 double edge_dest_p, edge_dest_q;
373 if (vertex_q[i] < 0) {
379 double edge_slope_p = edge_dest_p - edge_origin_p;
380 double edge_slope_q = edge_dest_q - edge_origin_q;
384 <<
" is inside and adjacent to a vertex outside the module.";
385 ldmx_log(trace) <<
"Working on edge with slope (" << edge_slope_p
386 <<
"," << edge_slope_q <<
")" <<
" and origin ("
387 << edge_origin_p <<
"," << edge_origin_q <<
")";
391 double projection_factor =
392 ((vertex_p[i] - edge_origin_p) * edge_slope_p +
393 (vertex_q[i] - edge_origin_q) * edge_slope_q) /
394 (edge_slope_p * edge_slope_p + edge_slope_q * edge_slope_q);
396 double proj_p = edge_origin_p + projection_factor * edge_slope_p;
397 double proj_q = edge_origin_q + projection_factor * edge_slope_q;
399 if (not isinside[up]) {
401 actual_p[num_vertices] = vertex_p[i];
402 actual_q[num_vertices] = vertex_q[i];
403 actual_p[num_vertices + 1] = proj_p;
404 actual_q[num_vertices + 1] = proj_q;
407 actual_p[num_vertices] = proj_p;
408 actual_q[num_vertices] = proj_q;
409 actual_p[num_vertices + 1] = vertex_p[i];
410 actual_q[num_vertices + 1] = vertex_q[i];
414 ldmx_log(trace) <<
"New Vertex " << i <<
" : (" << vertex_p[i]
415 <<
"," << vertex_q[i] <<
")";
418 actual_p[num_vertices] = vertex_p[i];
419 actual_q[num_vertices] = vertex_q[i];
426 for (
int i = 0; i < 6; i++) {
427 actual_p[i] = vertex_p[i];
428 actual_q[i] = vertex_q[i];
441 double p = (poly_bin->GetXMax() + poly_bin->GetXMin()) / 2.;
442 double q = (poly_bin->GetYMax() + poly_bin->GetYMin()) / 2.;
444 ldmx_log(trace) <<
" Copying poly with ID " << poly_bin->GetBinNumber()
445 <<
" and (p,q) (" << std::setprecision(2) << p <<
"," << q
456 ldmx_log(trace) <<
"Building cellModule position map";
460 double cell_x{cell_pq.first}, cell_y{cell_pq.second};
467 auto cell_rel_to_layer =
468 std::make_pair(module_xy.first + cell_x, module_xy.second + cell_y);
482 std::make_tuple(rel_to_layer.first + std::get<0>(layer_xyz),
483 rel_to_layer.second + std::get<1>(layer_xyz),
484 std::get<2>(layer_xyz));
507 ldmx_log(trace) <<
"Building Nearest and Next-Nearest Neighbor maps";
514 double dist = distance(probe_xyz, center_xyz);
516 nn_map_[center_id].push_back(probe_id);
518 nnn_map_[center_id].push_back(probe_id);
531 double x = fabs(cell_location.first);
532 double y = fabs(cell_location.second);
533 double r = sqrt(x * x + y * y);
534 double theta = (r > 1E-3) ? fabs(std::atan(y / x)) : 0;
539 sqrt(3.) *
module_r_max_ / (std::sin(theta) + sqrt(3.) * std::cos(theta));
544 ldmx_log(trace) << std::fixed << std::setprecision(2)
545 <<
"Checking if normXY=(" << normX <<
"," << normY
547 normX = fabs(normX), normY = fabs(normY);
548 double xvec = -1, yvec = -1. / sqrt(3);
549 double xref = 0.5, yref = sqrt(3) / 2.;
550 if ((normX > 1.) || (normY > yref)) {
551 ldmx_log(trace) <<
"They are outside quadrant.";
554 double dot_prod = (xvec * (normX - xref) + yvec * (normY - yref));
555 ldmx_log(trace) << std::fixed << std::setprecision(2)
556 <<
"They are inside quadrant. Dot product (>0 is inside): "
558 return (dot_prod > 0.);
Class that translates raw positions of ECal module hits into cells in a hexagonal readout.
Class encapsulating parameters for configuring a processor.
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Translation between real-space positions and cell IDs within the ECal.
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 > > nnn_map_
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.
std::tuple< double, double, double > getPosition(EcalID id) const
Get a cell's position from its ID number.
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.
bool layer_shift_odd_
shift odd layers
double module_r_max_
Center-to-Corner Radius of module hexagon [mm].
double cell_r_min_
Center-to-Flat Radius of cell hexagon [mm].
std::map< int, std::tuple< double, double, double > > layer_pos_xy_
Position of layer centers in world coordinates (uses layer ID as key)
double n_cell_r_height_
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::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]
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.
std::map< EcalID, std::vector< EcalID > > nn_map_
Map of cell ID to neighboring cells.
std::vector< double > layer_z_positions_
The layer Z postions are with respect to the front of the ECal [mm].
bool isInside(double normX, double normY) const
Determines if point (x,y), already normed to max hexagon radius, lies within a hexagon.
double si_thickness_
Thickness of the Si sensitive layer [mm].
bool layer_shift_odd_bilayer_
shift odd bi layers
void buildLayerMap()
Constructs the positions of the layers in world coordinates.
EcalGeometry(const framework::config::Parameters &ps)
Class constructor, for use only by the provider.
double cell_r_max_
Center-to-Corner Radius of cell hexagon [mm].
double ecal_front_z_
Front of ECal relative to world geometry [mm].
double gap_
Gap between module flat sides [mm].
double module_r_min_
Center-to-Flat Radius of module hexagon [mm].
TH2Poly cell_id_in_module_
Honeycomb Binning from ROOT.
void buildNeighborMaps()
Construts NNMap and NNNMap.
bool corners_side_up_
indicator of geometry orientation if true, flower shape's corners side (ie: side with two modules) is...
Extension of DetectorID providing access to ECal layers and cell numbers in a hex grid.
All classes in the ldmx-sw project use this namespace.