CMS 3D CMS Logo

TotemT2Segmentation.cc
Go to the documentation of this file.
1 
2 /****************************************************************************
3  *
4  * This is a part of TOTEM offline software.
5  * Author:
6  * Laurent Forthomme
7  * Arkadiusz Cwikla
8  *
9  ****************************************************************************/
10 
13 
14 #include "TH2D.h"
15 
16 TotemT2Segmentation::TotemT2Segmentation(size_t nbinsx, size_t nbinsy) : nbinsx_(nbinsx), nbinsy_(nbinsy) {
17  for (unsigned short arm = 0; arm <= CTPPSDetId::maxArm; ++arm)
18  for (unsigned short plane = 0; plane <= TotemT2DetId::maxPlane; ++plane)
19  for (unsigned short id = 0; id <= TotemT2DetId::maxChannel; ++id) {
20  const TotemT2DetId detid(arm, plane, id);
21  bins_map_[detid] = computeBins(detid);
22  }
23 }
24 
25 void TotemT2Segmentation::fill(TH2D* hist, const TotemT2DetId& detid, double value) {
26  if (bins_map_.count(detid) == 0)
27  throw cms::Exception("TotemT2Segmentation") << "Failed to retrieve list of bins for TotemT2DetId " << detid << ".";
28  if (static_cast<size_t>(hist->GetXaxis()->GetNbins()) != nbinsx_ ||
29  static_cast<size_t>(hist->GetYaxis()->GetNbins()) != nbinsy_)
30  throw cms::Exception("TotemT2Segmentation")
31  << "Trying to fill a summary plot with invalid number of bins. "
32  << "Should be of dimension (" << nbinsx_ << ", " << nbinsy_ << "), but has dimension ("
33  << hist->GetXaxis()->GetNbins() << ", " << hist->GetYaxis()->GetNbins() << ").";
34  for (const auto& bin : bins_map_.at(detid))
35  hist->Fill(bin.first, bin.second, value);
36 }
37 
38 std::vector<std::pair<short, short> > TotemT2Segmentation::computeBins(const TotemT2DetId& detid) const {
39  std::vector<std::pair<short, short> > bins;
40  // find the histogram centre
41  const auto ox = floor(nbinsx_ * 0.5), oy = floor(nbinsy_ * 0.5);
42  // compute the ellipse parameters
43  const auto ax = ceil(nbinsx_ * 0.5), by = ceil(nbinsy_ * 0.5);
44 
45  const float max_half_angle_rad = 0.3;
46 
47  // Without use of the Totem Geometry, follow an octant-division of channels,
48  // with even and odd planes alternating around the circle. Starting at 9AM
49  // and going clockwise around in increments of 1h30min=45deg, we have
50  // Ch0 on even planes, Ch0+odd, Ch1+even, Ch1+odd, up to Ch3+odd at 7:30PM
51 
52  const float tile_angle_rad = (180 - 45. / 2 - (detid.plane() % 2 ? 45 : 0) - (detid.channel() * 90)) * M_PI / 180.;
53  // Geometric way of associating a DetId to a vector<ix, iy> of bins given the size (nx_, ny_) of
54  // the TH2D('s) to be filled
55  for (size_t ix = 0; ix < nbinsx_; ++ix)
56  for (size_t iy = 0; iy < nbinsy_; ++iy) {
57  const auto ell_rad_norm = std::pow((ix - ox) / ax, 2) + std::pow((iy - oy) / by, 2);
58  if (ell_rad_norm < 1. && ell_rad_norm >= 0.1 &&
59  fabs(std::atan2(iy - oy, ix - ox) - tile_angle_rad) < max_half_angle_rad)
60  bins.emplace_back(ix, iy);
61  }
62 
63  return bins;
64 }
constexpr int32_t ceil(float num)
TotemT2Segmentation(size_t, size_t)
Detector ID class for Totem T2 detectors. Bits [19:31] : Base CTPPSDetId class attributes Bits [16:18...
Definition: TotemT2DetId.h:25
void fill(TH2D *, const TotemT2DetId &, double value=1.)
static constexpr uint32_t maxChannel
Definition: TotemT2DetId.h:35
constexpr int pow(int x)
Definition: conifer.h:24
static constexpr uint32_t maxPlane
Definition: TotemT2DetId.h:34
uint32_t plane() const
Definition: TotemT2DetId.h:44
static const uint32_t maxArm
Definition: CTPPSDetId.h:51
Definition: value.py:1
#define M_PI
std::unordered_map< TotemT2DetId, std::vector< std::pair< short, short > > > bins_map_
std::vector< std::pair< short, short > > computeBins(const TotemT2DetId &detid) const
uint32_t channel() const
Definition: TotemT2DetId.h:51