LDMX Software
EcalGeometry.cxx
2
3namespace ldmx {
4
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));
9}
10
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)));
20}
21
25static void rotate(double& p, double& q) {
26 double tmp{p};
27 p = -q;
28 q = tmp;
29}
30
34static void unrotate(double& p, double& q) {
35 double tmp{p};
36 p = q;
37 q = -tmp;
38}
39
41 : framework::ConditionsObject(EcalGeometry::CONDITIONS_OBJECT_NAME) {
42 layer_z_positions_ = ps.get<std::vector<double>>("layer_z_positions");
43 ecal_front_z_ = ps.get<double>("ecal_front_z");
44 module_r_min_ = ps.get<double>("module_min_r");
45 gap_ = ps.get<double>("gap");
46 n_cell_r_height_ = ps.get<double>("n_cell_r_height");
47 corners_side_up_ = ps.get<bool>("corners_side_up");
48 layer_shift_x_ = ps.get<double>("layer_shift_x");
49 layer_shift_y_ = ps.get<double>("layer_shift_y");
50 layer_shift_odd_ = ps.get<bool>("layer_shift_odd");
51 layer_shift_odd_bilayer_ = ps.get<bool>("layer_shift_odd_bilayer");
52 si_thickness_ = ps.get<double>("si_thickness");
53
55 EXCEPTION_RAISE("BadConf",
56 "Cannot shift both odd sensitive layers and odd bilayers");
57 }
58
59 module_r_max_ = module_r_min_ * (2 / sqrt(3));
61 cell_r_min_ = (sqrt(3.) / 2.) * cell_r_max_;
62
63 ldmx_log(debug) << "Building module map with gap " << std::setprecision(2)
64 << gap_ << ", nCellRHeight " << n_cell_r_height_
65 << ", min/max radii of cell " << cell_r_min_ << " / "
66 << cell_r_max_ << ", and module " << module_r_min_ << " / "
68
74
75 ldmx_log(trace) << "Geometry fully constructed";
76}
77
78EcalID EcalGeometry::getID(double x, double y, double z, bool fallible) const {
79 int layer_id{-1};
80 for (const auto& [lid, layer_xyz] : layer_pos_xy_) {
81 // check if the z coordinate is within the layer thickness
82 if (std::abs(std::get<2>(layer_xyz) - z) < si_thickness_) {
83 layer_id = lid;
84 break;
85 }
86 }
87 if (layer_id < 0) {
88 if (fallible) {
89 return ldmx::EcalID(0);
90 } else {
91 EXCEPTION_RAISE("BadConf", "z = " + std::to_string(z) +
92 " mm is not within any"
93 " of the configured layers.");
94 }
95 }
96 return getID(x, y, layer_id, fallible);
97}
98
99EcalID EcalGeometry::getID(double x, double y, int layer_id,
100 bool fallible) const {
101 // now assume we know the layer
102 // shift to center of layer
103 // and convert to flower coordinates
104 double p{x - std::get<0>(layer_pos_xy_.at(layer_id))},
105 q{y - std::get<1>(layer_pos_xy_.at(layer_id))};
106
107 // deduce module ID
108 // there are only 7 modules so we just loop through them
109 // all and pick out the module ID that we are inside of
110
111 int module_id{-1};
112 for (auto const& [mid, module_xy] : module_pos_xy_) {
113 double probe_x{p - module_xy.first}, probe_y{q - module_xy.second};
114 if (corners_side_up_) rotate(probe_x, probe_y);
115 if (isInside(probe_x / module_r_max_, probe_y / module_r_max_)) {
116 module_id = mid;
117 break;
118 }
119 }
120
121 if (module_id < 0) {
122 if (fallible) {
123 return ldmx::EcalID(0);
124 } else {
125 EXCEPTION_RAISE(
126 "BadConf",
127 TString::Format(
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)
132 .Data());
133 }
134 }
135
136 return getID(x, y, layer_id, module_id);
137}
138
139EcalID EcalGeometry::getID(double x, double y, int layer_id, int module_id,
140 bool fallible) const {
141 // now assume we know the layer and module
142 // shift to center of layer and then center of module
143 double p{x - std::get<0>(layer_pos_xy_.at(layer_id)) -
144 module_pos_xy_.at(module_id).first},
145 q{y - std::get<1>(layer_pos_xy_.at(layer_id)) -
146 module_pos_xy_.at(module_id).second};
147
148 // need to rotate
149 if (corners_side_up_) rotate(p, q);
150
151 // deduce cell ID
152 int cell_id = cell_id_in_module_.FindBin(p, q) - 1;
153
154 if (cell_id < 0) {
155 if (fallible) {
156 return ldmx::EcalID(0);
157 } else {
158 EXCEPTION_RAISE(
159 "BadConf",
160 TString::Format(
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)
165 .Data());
166 }
167 }
168
169 return EcalID(layer_id, module_id, cell_id);
170}
171
172std::tuple<double, double, double> EcalGeometry::getPosition(EcalID id) const {
173 return cell_global_pos_.at(id);
174}
175
176std::pair<double, double> EcalGeometry::getPositionInModule(int cell_id) const {
177 auto pq = cell_pos_in_module_.at(cell_id);
178
179 // going from (p,q) to (x,y) is a unrotate
180 if (corners_side_up_) unrotate(pq.first, pq.second);
181
182 return pq;
183}
184
186 ldmx_log(debug) << "Building layer map with " << layer_z_positions_.size()
187 << " layers";
189 if (layer_shift_odd_) {
190 ldmx_log(debug) << "Shifting odd layers by (x,y) = (" << layer_shift_x_
191 << ", " << layer_shift_y_ << ") mm";
192 } else {
193 ldmx_log(debug) << "Shifting odd bilayers by (x,y) = (" << layer_shift_x_
194 << ", " << layer_shift_y_ << ") mm";
195 }
196 }
197
198 for (std::size_t i_layer{0}; i_layer < layer_z_positions_.size(); ++i_layer) {
199 // default is centered on z-axis
200 double x{0}, y{0}, z{ecal_front_z_ + layer_z_positions_.at(i_layer)};
201 if (layer_shift_odd_ and (i_layer % 2 == 1)) {
202 x += layer_shift_x_;
203 y += layer_shift_y_;
204 } else if (layer_shift_odd_bilayer_ and ((i_layer / 2) % 2 == 1)) {
205 x += layer_shift_x_;
206 y += layer_shift_y_;
207 }
208 ldmx_log(trace) << " Layer " << i_layer << " has center at (" << x
209 << ", " << y << ", " << z << ") mm";
210
211 layer_pos_xy_[i_layer] = std::make_tuple(x, y, z);
212 }
213}
214
216 static const double c_pi = 3.14159265358979323846;
217
218 // the center module (module_id == 0) has always been (and will always be?)
219 // centered with respect to the layer position
220 module_pos_xy_[0] = std::pair<double, double>(0., 0.);
221
222 // for flat-side-up designs (v12 and earlier), the modules are numbered 1 on
223 // positive y-axis and then counter-clockwise until 6
224 //
225 // for corner-side-up designes (v13 and later), the modules are numbered 1 on
226 // positive x-axis and then counter-clockwise until 6.
227 for (unsigned id = 1; id < 7; id++) {
228 // flat-side-up
229 double x = (2. * module_r_min_ + gap_) * sin((id - 1) * (c_pi / 3.));
230 double y = (2. * module_r_min_ + gap_) * cos((id - 1) * (c_pi / 3.));
231 if (corners_side_up_) {
232 // re-calculating to make sure centers match GDML
233 x = (2. * module_r_min_ + gap_) * cos((id - 1) * (c_pi / 3.));
234 y = -(2. * module_r_min_ + gap_) * sin((id - 1) * (c_pi / 3.));
235 }
236 module_pos_xy_[id] = std::pair<double, double>(x, y);
237 ldmx_log(trace) << " Module " << id << " is centered at (x,y) = " << "("
238 << x << ", " << y << ") mm";
239 }
240}
241
262 TH2Poly grid_map;
263
264 // make hexagonal grid [boundary is rectangle] larger than the module
265 double grid_min_p = 0., grid_min_q = 0.; // start at the origin
266 int num_p_cells = 0, num_q_cells = 0;
267
268 // first x-cell is only a half
269 grid_min_p -= cell_r_min_;
270 num_p_cells++;
271 while (grid_min_p > -1 * module_r_max_) {
272 // decrement x by cell center-to-flat diameter
273 grid_min_p -= 2 * cell_r_min_;
274 num_p_cells++;
275 }
276 while (grid_min_q > -1 * module_r_min_) {
277 // decrement y by cell center-to-corner radius
278 // alternate between a full corner-to-corner diameter
279 // and a side of a cell (center-to-corner radius)
280 if (num_q_cells % 2 == 0)
281 grid_min_q -= 1 * cell_r_max_;
282 else
283 grid_min_q -= 2 * cell_r_max_;
284 num_q_cells++;
285 }
286 // only counted one half of the cells
287 num_p_cells *= 2;
288 num_q_cells *= 2;
289
290 grid_map.Honeycomb(grid_min_p, grid_min_q, cell_r_max_, num_p_cells,
291 num_q_cells);
292
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 << ")";
299
300 // copy cells lying within module boundaries to a module grid
301 TListIter next(grid_map.GetBins()); // a TH2Poly is a TList of TH2PolyBin
302 TH2PolyBin* poly_bin = 0;
303 TGraph* poly = 0; // a polygon returned by TH2Poly is a TGraph
304 // cells_in_module_ IDs go from 0 to N-1, not equal to original grid cell ID
305 int cell_id = 0;
306 while ((poly_bin = (TH2PolyBin*)next())) {
307 // these bins are coming from the honeycomb
308 // grid so we assume that they are regular
309 // hexagons i.e. has 6 vertices
310 poly = (TGraph*)poly_bin->GetPolygon();
311
312 // decide whether to copy polygon to new map.
313 // use all vertices in case of cut-off edge polygons.
314 int num_vertices_inside = 0;
315 double vertex_p[6], vertex_q[6]; // vertices of the cell
316 bool isinside[6]; // which vertices are inside
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];
322
323 isinside[i] =
324 isInside(vertex_p[i] / module_r_max_, vertex_q[i] / module_r_max_);
325 if (isinside[i]) num_vertices_inside++;
326 }
327
328 if (num_vertices_inside > 1) {
329 // Include this cell if more than one of its vertices is inside the module
330 // hexagon
331 double actual_p[8], actual_q[8];
332 int num_vertices{0};
333 if (num_vertices_inside < 6) {
334 // This cell is stradling the edge of the module
335 // and is NOT cleanly cut by module edge
336
337 ldmx_log(trace) << " Polygon " << cell_id
338 << " has vertices poking out of module hexagon.";
339
340 // loop through vertices
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])) {
345 // this vertex is inside the module hexagon and is adjacent to a
346 // vertex outside
347 // ==> project this vertex onto the nearest edge of the module
348 // hexagon
349
350 // determine which side of hexagon we should project onto
351 double edge_origin_p, edge_origin_q;
352 double edge_dest_p, edge_dest_q;
353 if (vertex_p[i] < -module_r_max_ / 2.) {
354 // sloped edge on negative-x side
355 edge_origin_p = -1. * module_r_max_;
356 edge_origin_q = 0.;
357 edge_dest_p = -0.5 * module_r_max_;
358 edge_dest_q = module_r_min_;
359 } else if (vertex_p[i] > module_r_max_ / 2.) {
360 // sloped edge on positive-x side
361 edge_origin_p = 0.5 * module_r_max_;
362 edge_origin_q = module_r_min_;
363 edge_dest_p = module_r_max_;
364 edge_dest_q = 0.;
365 } else {
366 // flat edge at top
367 edge_origin_p = 0.5 * module_r_max_;
368 edge_origin_q = module_r_min_;
369 edge_dest_p = -0.5 * module_r_max_;
370 edge_dest_q = module_r_min_;
371 }
372 // flip to bottom half if below x-axis
373 if (vertex_q[i] < 0) {
374 edge_dest_q *= -1;
375 edge_origin_q *= -1;
376 }
377
378 // get edge slope vector
379 double edge_slope_p = edge_dest_p - edge_origin_p;
380 double edge_slope_q = edge_dest_q - edge_origin_q;
381
382 ldmx_log(trace)
383 << "Vertex " << i
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 << ")";
388
389 // project vertices adjacent to the vertex outside the module onto
390 // the module edge
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);
395
396 double proj_p = edge_origin_p + projection_factor * edge_slope_p;
397 double proj_q = edge_origin_q + projection_factor * edge_slope_q;
398
399 if (not isinside[up]) {
400 // the next point is outside
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;
405 } else {
406 // the previous point was outside
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];
411 }
412 num_vertices += 2;
413
414 ldmx_log(trace) << "New Vertex " << i << " : (" << vertex_p[i]
415 << "," << vertex_q[i] << ")";
416
417 } else {
418 actual_p[num_vertices] = vertex_p[i];
419 actual_q[num_vertices] = vertex_q[i];
420 num_vertices++;
421 } // should we project or not
422 } // loop through vertices
423 } else {
424 // all 6 inside, just copy the vertices over
425 num_vertices = 6;
426 for (int i = 0; i < 6; i++) {
427 actual_p[i] = vertex_p[i];
428 actual_q[i] = vertex_q[i];
429 }
430 } // if numVerticesInside is less than 5
431
432 // TH2Poly needs to have its own copy of the polygon TGraph
433 // otherwise, we get a seg fault when EcalGeometry is destructed
434 // because the polygon that was copied over from gridMap is deleted at
435 // the end of this function
436 cell_id_in_module_.AddBin(num_vertices, actual_p, actual_q);
437
441 double p = (poly_bin->GetXMax() + poly_bin->GetXMin()) / 2.;
442 double q = (poly_bin->GetYMax() + poly_bin->GetYMin()) / 2.;
443
444 ldmx_log(trace) << " Copying poly with ID " << poly_bin->GetBinNumber()
445 << " and (p,q) (" << std::setprecision(2) << p << "," << q
446 << ")";
447 // save cell location as center of ENTIRE hexagon
448 cell_pos_in_module_[cell_id] = std::pair<double, double>(p, q);
449 ++cell_id; // incrememnt cell ID
450 } // if num vertices inside is > 1
451 } // loop over larger grid spanning module hexagon
452 return;
453}
454
456 ldmx_log(trace) << "Building cellModule position map";
458 for (auto const& [module_id, module_xy] : module_pos_xy_) {
459 for (auto const& [cell_id, cell_pq] : cell_pos_in_module_) {
460 double cell_x{cell_pq.first}, cell_y{cell_pq.second};
461 // convert from (p,q) to (x,y) space
462 // when the corners are not up, x = p and y = q
463 // so no transformation needs to be done
464 if (corners_side_up_) unrotate(cell_x, cell_y);
465
466 // calculate cell's pq relative to entire layer center
467 auto cell_rel_to_layer =
468 std::make_pair(module_xy.first + cell_x, module_xy.second + cell_y);
469
470 // now add the layer-center values to get the global position
471 // of the cell
472 cell_pos_in_layer_[EcalID(0, module_id, cell_id)] = cell_rel_to_layer;
473 }
474 }
475
477 for (auto const& [layer_id, layer_xyz] : layer_pos_xy_) {
478 for (auto const& [flat_id, rel_to_layer] : cell_pos_in_layer_) {
479 // now add the layer-center values to get the global position
480 // of the cell
481 cell_global_pos_[EcalID(layer_id, flat_id.module(), flat_id.cell())] =
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));
485 }
486 }
487
488 ldmx_log(trace) << "Cell module map contains " << cell_global_pos_.size()
489 << " entries. ";
490 return;
491}
492
507 ldmx_log(trace) << "Building Nearest and Next-Nearest Neighbor maps";
508
509 nn_map_.clear();
510 nnn_map_.clear();
511 for (auto const& [center_id, center_xyz] : cell_pos_in_layer_) {
512 for (auto const& [probe_id, probe_xyz] : cell_pos_in_layer_) {
514 double dist = distance(probe_xyz, center_xyz);
515 if (dist > 1 * cell_r_min_ && dist <= 3. * cell_r_min_) {
516 nn_map_[center_id].push_back(probe_id);
517 } else if (dist > 3. * cell_r_min_ && dist <= 4.5 * cell_r_min_) {
518 nnn_map_[center_id].push_back(probe_id);
519 }
520 }
521 // Keeping this commented out until the SimIDs string method is implemented
522 // ldmx_log(debug) << " Found " << nn_map_[center_id].size() << " NN and "
523 // << nnn_map_[center_id].size() << " NNN for cell " << center_id;
524 }
525 return;
526}
527
529 // https://math.stackexchange.com/questions/1210572/find-the-distance-to-the-edge-of-a-hexagon
530 std::pair<double, double> cell_location = cell_pos_in_module_.at(id.cell());
531 double x = fabs(cell_location.first); // bring to first quadrant
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;
535 if (x < module_r_max_ / 2.)
536 return (module_r_min_ -
537 y); // closest line is straight vertical to top edge
538 double dist =
539 sqrt(3.) * module_r_max_ / (std::sin(theta) + sqrt(3.) * std::cos(theta));
540 return dist;
541}
542
543bool EcalGeometry::isInside(double normX, double normY) const {
544 ldmx_log(trace) << std::fixed << std::setprecision(2)
545 << "Checking if normXY=(" << normX << "," << normY
546 << ") is inside.";
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.";
552 return false;
553 }
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): "
557 << dot_prod;
558 return (dot_prod > 0.);
559}
560
561} // namespace ldmx
Class that translates raw positions of ECal module hits into cells in a hexagonal readout.
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
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.
Definition EcalID.h:20
All classes in the ldmx-sw project use this namespace.