ATLAS Offline Software
Loading...
Searching...
No Matches
CaloGeometryLookup.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 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
14//#include "ISF_FastCaloSimEvent/ICaloGeometry.h"
17
19class TCanvas;
20class TGraphErrors;
21
22typedef std::map< Identifier , const CaloDetDescrElement* > t_cellmap;
23typedef std::map< double , const CaloDetDescrElement* > t_eta_cellmap;
24
26 public:
27 CaloGeometryLookup(int ind=0);
29
30 bool IsCompatible(const CaloDetDescrElement* cell);
31 void add(const CaloDetDescrElement* cell);
32 t_cellmap::size_type size() const {return m_cells.size();};
33 int index() const {return m_index;};
34 void set_index(int ind) {m_index=ind;};
35 void post_process();
38 //void CalculateTransformation();
39
40 float mineta() const {return m_mineta;};
41 float maxeta() const {return m_maxeta;};
42 float minphi() const {return m_minphi;};
43 float maxphi() const {return m_maxphi;};
44
45 float mineta_raw() const {return m_mineta_raw;};
46 float maxeta_raw() const {return m_maxeta_raw;};
47 float minphi_raw() const {return m_minphi_raw;};
48 float maxphi_raw() const {return m_maxphi_raw;};
49
50#ifdef USE_GPU
51 //the geometry for GPU EDM are constructed with following functions
52 float mineta_correction() const {return m_mineta_correction;};
53 float maxeta_correction() const {return m_maxeta_correction;};
54 float minphi_correction() const {return m_minphi_correction;};
55 float maxphi_correction() const {return m_maxphi_correction;};
56#endif
57
58 float minx() const {return m_mineta;};
59 float maxx() const {return m_maxeta;};
60 float miny() const {return m_minphi;};
61 float maxy() const {return m_maxphi;};
62
63 float minx_raw() const {return m_mineta_raw;};
64 float maxx_raw() const {return m_maxeta_raw;};
65 float miny_raw() const {return m_minphi_raw;};
66 float maxy_raw() const {return m_maxphi_raw;};
67
68 const MeanAndRMS& deta() {return m_deta;};
69 const MeanAndRMS& dphi() {return m_dphi;};
70 float mindeta() {return m_mindeta;};
71 float mindphi() {return m_mindphi;};
72 const MeanAndRMS& dx() {return m_deta;};
73 const MeanAndRMS& dy() {return m_dphi;};
74 float mindx() {return m_mindeta;};
75 float mindy() {return m_mindphi;};
76
81
82 int cell_grid_eta() const {return m_cell_grid_eta;};
83 int cell_grid_phi() const {return m_cell_grid_phi;};
85
86 virtual const CaloDetDescrElement* getDDE(float eta,float phi,float* distance=0,int* steps=0);
87
88#ifdef USE_GPU
89 //used by FCS-GPU
90 float xy_grid_adjustment_factor() const {return m_xy_grid_adjustment_factor;};
91 float deta_double() const {return m_deta_double;};
92 float dphi_double() const {return m_dphi_double;};
93 const std::vector< std::vector< const CaloDetDescrElement* > > * cell_grid(){ return & m_cell_grid ; } ;
94#endif
95
96 protected:
97 float neta_double() {return (maxeta_raw()-mineta_raw())/deta().mean();};
98 float nphi_double() {return (maxphi_raw()-minphi_raw())/dphi().mean();};
99 Int_t neta() {return TMath::Nint( neta_double() );};
100 Int_t nphi() {return TMath::Nint( nphi_double() );};
101
102 //FCal is not sufficiently regular to use submaps with regular mapping
103 float nx_double() {return (maxx_raw()-minx_raw())/mindx();};
104 float ny_double() {return (maxy_raw()-miny_raw())/mindy();};
105 Int_t nx() {return TMath::Nint(TMath::Ceil( nx_double() ));};
106 Int_t ny() {return TMath::Nint(TMath::Ceil( ny_double() ));};
107
109
110 int raw_eta_position_to_index(float eta_raw) const {return TMath::Floor((eta_raw-mineta_raw())/m_deta_double);};
111 int raw_phi_position_to_index(float phi_raw) const {return TMath::Floor((phi_raw-minphi_raw())/m_dphi_double);};
112 bool index_range_adjust(int& ieta,int& iphi) const;
113 float calculate_distance_eta_phi(const CaloDetDescrElement* DDE,float eta,float phi,float& dist_eta0,float& dist_phi0) const;
114
117 std::vector< std::vector< const CaloDetDescrElement* > > m_cell_grid;
122
126};
127#endif
128
129
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)
double mean() const
Definition MeanAndRMS.h:22