ATLAS Offline Software
Loading...
Searching...
No Matches
CaloGeometryLookup.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef ISF_FASTCALOSIMPARAMETRIZATION_CALOGEOMETRYLOOKUP_H
6#define ISF_FASTCALOSIMPARAMETRIZATION_CALOGEOMETRYLOOKUP_H
7
8#include <TMath.h>
9
10#include <vector>
11#include <map>
12#include <iostream>
13
16
18class TCanvas;
19class TGraphErrors;
20
21typedef std::map< Identifier , const CaloDetDescrElement* > t_cellmap;
22typedef std::map< double , const CaloDetDescrElement* > t_eta_cellmap;
23
25 public:
26 CaloGeometryLookup(int ind=0);
28
29 bool IsCompatible(const CaloDetDescrElement* cell);
30 void add(const CaloDetDescrElement* cell);
31 t_cellmap::size_type size() const {return m_cells.size();};
32 int index() const {return m_index;};
33 void set_index(int ind) {m_index=ind;};
34 void post_process();
37 //void CalculateTransformation();
38
39 float mineta() const {return m_mineta;};
40 float maxeta() const {return m_maxeta;};
41 float minphi() const {return m_minphi;};
42 float maxphi() const {return m_maxphi;};
43
44 float mineta_raw() const {return m_mineta_raw;};
45 float maxeta_raw() const {return m_maxeta_raw;};
46 float minphi_raw() const {return m_minphi_raw;};
47 float maxphi_raw() const {return m_maxphi_raw;};
48
49#ifdef USE_GPU
50 //the geometry for GPU EDM are constructed with following functions
51 float mineta_correction() const {return m_mineta_correction;};
52 float maxeta_correction() const {return m_maxeta_correction;};
53 float minphi_correction() const {return m_minphi_correction;};
54 float maxphi_correction() const {return m_maxphi_correction;};
55#endif
56
57 float minx() const {return m_mineta;};
58 float maxx() const {return m_maxeta;};
59 float miny() const {return m_minphi;};
60 float maxy() const {return m_maxphi;};
61
62 float minx_raw() const {return m_mineta_raw;};
63 float maxx_raw() const {return m_maxeta_raw;};
64 float miny_raw() const {return m_minphi_raw;};
65 float maxy_raw() const {return m_maxphi_raw;};
66
67 const MeanAndRMS& deta() {return m_deta;};
68 const MeanAndRMS& dphi() {return m_dphi;};
69 float mindeta() {return m_mindeta;};
70 float mindphi() {return m_mindphi;};
71 const MeanAndRMS& dx() {return m_deta;};
72 const MeanAndRMS& dy() {return m_dphi;};
73 float mindx() {return m_mindeta;};
74 float mindy() {return m_mindphi;};
75
80
81 int cell_grid_eta() const {return m_cell_grid_eta;};
82 int cell_grid_phi() const {return m_cell_grid_phi;};
84
85 virtual const CaloDetDescrElement* getDDE(float eta,float phi,float* distance=0,int* steps=0);
86
87#ifdef USE_GPU
88 //used by FCS-GPU
89 float xy_grid_adjustment_factor() const {return m_xy_grid_adjustment_factor;};
90 float deta_double() const {return m_deta_double;};
91 float dphi_double() const {return m_dphi_double;};
92 const std::vector< std::vector< const CaloDetDescrElement* > > * cell_grid(){ return & m_cell_grid ; } ;
93#endif
94
95 protected:
96 float neta_double() {
97 if (auto m = deta().mean();m!=0.)[[likely]]{
98 return (maxeta_raw()-mineta_raw())/m;
99 }
100 return 0.f;
101 }
102
103 float nphi_double() {
104 if(auto m = dphi().mean(); m!=0.)[[likely]]{
105 return (maxphi_raw()-minphi_raw())/m;
106 }
107 return 0.f;
108 }
109
110 Int_t neta() {return TMath::Nint( neta_double() );};
111 Int_t nphi() {return TMath::Nint( nphi_double() );};
112
113 //FCal is not sufficiently regular to use submaps with regular mapping
114 float nx_double() {return (maxx_raw()-minx_raw())/mindx();};
115 float ny_double() {return (maxy_raw()-miny_raw())/mindy();};
116 Int_t nx() {return TMath::Nint(TMath::Ceil( nx_double() ));};
117 Int_t ny() {return TMath::Nint(TMath::Ceil( ny_double() ));};
118
120
121 int raw_eta_position_to_index(float eta_raw) const {return TMath::Floor((eta_raw-mineta_raw())/m_deta_double);};
122 int raw_phi_position_to_index(float phi_raw) const {return TMath::Floor((phi_raw-minphi_raw())/m_dphi_double);};
123 bool index_range_adjust(int& ieta,int& iphi) const;
124 float calculate_distance_eta_phi(const CaloDetDescrElement* DDE,float eta,float phi,float& dist_eta0,float& dist_phi0) const;
125
128 std::vector< std::vector< const CaloDetDescrElement* > > m_cell_grid;
133
137};
138#endif
139
140
const boost::regex ref(r_ef)
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
std::map< Identifier, const CaloDetDescrElement * > t_cellmap
std::map< double, const CaloDetDescrElement * > t_eta_cellmap
This class groups all DetDescr information related to a CaloCell.
const MeanAndRMS & deta()
std::vector< std::vector< const CaloDetDescrElement * > > m_cell_grid
const MeanAndRMS & dphi()
bool IsCompatible(const CaloDetDescrElement *cell)
int raw_phi_position_to_index(float phi_raw) const
t_cellmap::size_type size() const
void add(const CaloDetDescrElement *cell)
virtual const CaloDetDescrElement * getDDE(float eta, float phi, float *distance=0, int *steps=0)
bool index_range_adjust(int &ieta, int &iphi) const
bool has_overlap(CaloGeometryLookup *ref)
const MeanAndRMS & x_correction()
const MeanAndRMS & dy()
const MeanAndRMS & eta_correction()
const MeanAndRMS & y_correction()
const MeanAndRMS & phi_correction()
int raw_eta_position_to_index(float eta_raw) const
float calculate_distance_eta_phi(const CaloDetDescrElement *DDE, float eta, float phi, float &dist_eta0, float &dist_phi0) const
virtual ~CaloGeometryLookup()
const MeanAndRMS & dx()
void merge_into_ref(CaloGeometryLookup *ref)
void set_xy_grid_adjustment_factor(float factor)
void mean(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
#define likely(x)