ATLAS Offline Software
MuonSectorMapping.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #ifndef MUON_MUONSECTORMAPPING_H
6 #define MUON_MUONSECTORMAPPING_H
7 
8 #include <math.h>
9 
10 #include <array>
11 #include <iostream>
12 #include <vector>
13 
15 
17 
18 namespace Muon {
19 
21  public:
23  bool insideSector(int sector, double phi) const;
24 
26  int getSector(double phi) const;
27 
29  void getSectors(double phi, std::vector<int>& sectors) const;
30 
32  double sectorPhi(int sector) const;
33  /*** Returns whether the sector is a small sector */
34  static bool isSmall(int sector) {
35  return sector % 2 == 0;
36  }
37 
39  double transformPhiToSector(double phi, int sector, bool toSector = true) const;
40 
42  double transformRToSector(double r, double phi, int sector, bool toSector = true) const;
43 
45  double sectorOverlapPhi(int sector1, int sector2) const;
46 
48  bool closeToSectorBoundary(double phi) const;
49 
51  double transformRToNeighboringSector(double r, int sectorHit, int sectorTarget) const;
52 
54  double sectorSize(int sector) const;
55 
57  double sectorWidth(int sector) const;
58 
59  private:
60  static constexpr double s_oneEightsOfPi{M_PI / 8.}; // pi/8
61  static constexpr double s_inverseOneEightsOfPi{8. / M_PI}; // 8/pi
62  static constexpr std::array<double, 2> s_sectorSize{0.4 * s_oneEightsOfPi, 0.6 * s_oneEightsOfPi}; // side of a sector in radiants
63  static constexpr double s_sectorOverlap{0.1 * s_oneEightsOfPi}; // size of the overlap between small and large sectors
64  static bool s_debug ATLAS_THREAD_SAFE; // if true addition cout is enabled
65  };
66 
67  inline double MuonSectorMapping::sectorSize(int sector) const {
68  const int idx = sector % 2;
69  return s_sectorSize[idx];
70  }
71 
72  inline double MuonSectorMapping::sectorWidth(int sector) const { return sectorSize(sector) + s_sectorOverlap; }
73 
74  inline bool MuonSectorMapping::insideSector(int sector, double phi) const {
75  double phiInSec = transformPhiToSector(phi, sector);
76  if (phiInSec < -sectorWidth(sector)) return false;
77  if (phiInSec > sectorWidth(sector)) return false;
78  return true;
79  }
80 
81  inline double MuonSectorMapping::sectorPhi(int sector) const {
82  if (sector < 10) return M_PI * (sector - 1) / 8.;
83  return -M_PI * (2 - (sector - 1) / 8.);
84  }
85 
86  inline int MuonSectorMapping::getSector(double phi) const {
87  // remap phi onto sector structure
88  double val = (phi + sectorSize(1)) * s_inverseOneEightsOfPi; // convert to value between -8 and 8, shift by width first sector
89  if (val < 0) val += 16;
90  int sliceIndex = static_cast<int>(val / 2); // large/small wedge
91  double valueInSlice = val - 2 * sliceIndex;
92  int sector = 2 * sliceIndex + 1;
93  if (valueInSlice > 1.2) ++sector;
94  return sector;
95  }
96 
97  inline void MuonSectorMapping::getSectors(double phi, std::vector<int>& sectors) const {
98  int sector = getSector(phi);
99  int sectorNext = sector != 16 ? sector + 1 : 1;
100  int sectorBefore = sector != 1 ? sector - 1 : 16;
101  if (insideSector(sectorBefore, phi)) sectors.push_back(sectorBefore);
102  if (insideSector(sector, phi)) sectors.push_back(sector);
103  if (insideSector(sectorNext, phi)) sectors.push_back(sectorNext);
104  }
105 
106  inline bool MuonSectorMapping::closeToSectorBoundary(double phi) const {
107  std::vector<int> sectors;
108  getSectors(phi, sectors);
109  return sectors.size() > 1;
110  }
111 
112  inline double MuonSectorMapping::transformPhiToSector(double phi, int sector, bool toSector) const {
113  double sign = toSector ? -1 : 1;
114  double dphi = phi + sign * sectorPhi(sector);
115  if (dphi > M_PI) dphi -= 2 * M_PI;
116  if (dphi < -M_PI) dphi += 2 * M_PI;
117  return dphi;
118  }
119 
120  inline double MuonSectorMapping::sectorOverlapPhi(int sector1, int sector2) const {
121  if (sector1 == sector2) return sectorPhi(sector1);
122 
123  int s1 = sector1 < sector2 ? sector1 : sector2;
124  int s2 = sector1 > sector2 ? sector1 : sector2;
125  if (s2 == 16 && s1 == 1) {
126  s1 = 16;
127  s2 = 1;
128  } else if (std::abs(s1 - s2) > 1) {
129  if (s_debug) std::cout << " bad sector combination: not neighbouring " << s1 << " " << s2 << std::endl;
130  return 0;
131  }
132 
133  double phi1 = sectorPhi(s1);
134  double phio1 = phi1 + sectorSize(s1);
135  if (phio1 > M_PI) phio1 -= 2 * M_PI;
136 
137  return phio1;
138  }
139 
140  inline double MuonSectorMapping::transformRToSector(double r, double phi, int sector, bool toSector) const {
141  double dphi = transformPhiToSector(phi, sector);
142  if (std::abs(dphi) > 0.3 && s_debug) {
143  std::cout << " large dphi detected!!: phi " << phi << " sector " << sector << " phi " << sectorPhi(sector) << " " << phi
144  << " dphi " << dphi << std::endl;
145  }
146  if (toSector)
147  return r * std::cos(dphi);
148  else
149  return r / std::cos(dphi);
150  }
151 
152  inline double MuonSectorMapping::transformRToNeighboringSector(double r, int sectorHit, int sector) const {
153  double phi = sectorPhi(sectorHit);
154  double phio = sectorOverlapPhi(sectorHit, sector);
155  double dphio = phi - phio;
156  if (dphio < -M_PI) dphio += 2 * M_PI;
157  double redge = r / std::cos(dphio);
158  double phi_sec = sectorPhi(sector);
159  double dphi = phio - phi_sec;
160  if (dphi < -M_PI) dphi += 2 * M_PI;
161  if (std::abs(dphi) > 0.3 && s_debug) {
162  std::cout << " large dphi detected!!: sector " << sector << " of hit " << sectorHit << " phi ref sector " << phi_sec << " hit "
163  << phi << " dphi " << dphi << std::endl;
164  }
165  return redge * std::cos(dphi);
166  }
167 
168 } // namespace Muon
169 
170 #endif
beamspotman.r
def r
Definition: beamspotman.py:676
Muon::MuonSectorMapping::sectorSize
double sectorSize(int sector) const
sector size (exclusive) in radians
Definition: MuonSectorMapping.h:67
Muon::MuonSectorMapping::sectorWidth
double sectorWidth(int sector) const
sector width (with overlap) in radians
Definition: MuonSectorMapping.h:72
Muon::MuonSectorMapping::transformRToSector
double transformRToSector(double r, double phi, int sector, bool toSector=true) const
expresses a radial position from and to the sector coordinate frame, the phi position should always b...
Definition: MuonSectorMapping.h:140
Muon::MuonSectorMapping::closeToSectorBoundary
bool closeToSectorBoundary(double phi) const
checks whether the phi position is close to a sector boundary
Definition: MuonSectorMapping.h:106
Muon::MuonSectorMapping::getSectors
void getSectors(double phi, std::vector< int > &sectors) const
returns the main sector plus neighboring if the phi position is in an overlap region
Definition: MuonSectorMapping.h:97
M_PI
#define M_PI
Definition: ActiveFraction.h:11
Muon::MuonSectorMapping::transformPhiToSector
double transformPhiToSector(double phi, int sector, bool toSector=true) const
transforms a phi position from and to the sector coordinate system in radians
Definition: MuonSectorMapping.h:112
drawFromPickle.cos
cos
Definition: drawFromPickle.py:36
Muon
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
Definition: TrackSystemController.h:45
Muon::MuonSectorMapping::s_sectorSize
static constexpr std::array< double, 2 > s_sectorSize
Definition: MuonSectorMapping.h:62
Muon::MuonSectorMapping::getSector
int getSector(double phi) const
returns the sector corresponding to the phi position
Definition: MuonSectorMapping.h:86
ATLAS_CHECK_FILE_THREAD_SAFETY
ATLAS_CHECK_FILE_THREAD_SAFETY
Definition: MuonSectorMapping.h:16
Muon::MuonSectorMapping::insideSector
bool insideSector(int sector, double phi) const
checks whether the phi position is consistent with sector
Definition: MuonSectorMapping.h:74
Muon::MuonSectorMapping::s_oneEightsOfPi
static constexpr double s_oneEightsOfPi
Definition: MuonSectorMapping.h:60
sign
int sign(int a)
Definition: TRT_StrawNeighbourSvc.h:107
Muon::MuonSectorMapping::sectorPhi
double sectorPhi(int sector) const
returns the centeral phi position of a sector in radians
Definition: MuonSectorMapping.h:81
Muon::MuonSectorMapping::isSmall
static bool isSmall(int sector)
Definition: MuonSectorMapping.h:34
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
ReadCellNoiseFromCoolCompare.s2
s2
Definition: ReadCellNoiseFromCoolCompare.py:379
LArNewCalib_DelayDump_OFC_Cali.idx
idx
Definition: LArNewCalib_DelayDump_OFC_Cali.py:69
Muon::MuonSectorMapping::s_inverseOneEightsOfPi
static constexpr double s_inverseOneEightsOfPi
Definition: MuonSectorMapping.h:61
Muon::MuonSectorMapping
Definition: MuonSectorMapping.h:20
Muon::MuonSectorMapping::transformRToNeighboringSector
double transformRToNeighboringSector(double r, int sectorHit, int sectorTarget) const
transform a radial position from one sector frame into another
Definition: MuonSectorMapping.h:152
Muon::MuonSectorMapping::s_sectorOverlap
static constexpr double s_sectorOverlap
Definition: MuonSectorMapping.h:63
Muon::MuonSectorMapping::sectorOverlapPhi
double sectorOverlapPhi(int sector1, int sector2) const
returns the phi position of the overlap between the two sectors (which have to be neighboring) in rad...
Definition: MuonSectorMapping.h:120
checker_macros.h
Define macros for attributes used to control the static checker.
Muon::MuonSectorMapping::ATLAS_THREAD_SAFE
static bool s_debug ATLAS_THREAD_SAFE
Definition: MuonSectorMapping.h:64