Project
Loading...
Searching...
No Matches
Geometry.h
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
12#ifndef ALICEO2_EMCAL_GEOMETRY_H_
13#define ALICEO2_EMCAL_GEOMETRY_H_
14
15#include <array>
16#include <string>
17#include <string_view>
18#include <tuple>
19#include <vector>
20
21#include <RStringView.h>
22#include <RtypesCore.h>
23#include <TGeoMatrix.h>
24#include <TMath.h>
25#include <TParticle.h>
26#include <TVector3.h>
27
30#include "MathUtils/Cartesian.h"
31
32namespace o2
33{
34namespace emcal
35{
36class ShishKebabTrd1Module;
37
42{
43 public:
47 Geometry() = default;
48
62 explicit Geometry(const std::string_view name, const std::string_view mcname = "", const std::string_view mctitle = "");
63
65 Geometry(const Geometry& geom);
66
68 ~Geometry();
69
71 Geometry& operator=(const Geometry& rvalue);
72
76 static Geometry* GetInstance();
77
87 static Geometry* GetInstance(const std::string_view name, const std::string_view mcname = "TGeant3",
88 const std::string_view mctitle = "");
89
96 static Geometry* GetInstanceFromRunNumber(Int_t runNumber, const std::string_view = "",
97 const std::string_view mcname = "TGeant3",
98 const std::string_view mctitle = "");
99
105 void DefineSamplingFraction(const std::string_view mcname = "", const std::string_view mctitle = "");
106
108 // General
109 //
110
111 const std::string& GetName() const { return mGeoName; }
112
113 static const std::string& GetDefaultGeometryName() { return DEFAULT_GEOMETRY; }
114
115 static Bool_t IsInitialized() { return Geometry::sGeom != nullptr; }
116
123
124 const std::vector<ShishKebabTrd1Module>& GetShishKebabTrd1Modules() const { return mShishKebabTrd1Modules; }
125
128 const ShishKebabTrd1Module& GetShishKebabModule(Int_t neta) const;
129
135 Bool_t Impact(const TParticle* particle) const;
136
147 void ImpactOnEmcal(const math_utils::Point3D<double>& vtx, Double_t theta, Double_t phi, Int_t& absId, math_utils::Point3D<double>& vimpact) const;
148
154 Bool_t IsInEMCAL(const math_utils::Point3D<double>& pnt) const;
155
161 Bool_t IsInDCAL(const math_utils::Point3D<double>& pnt) const;
162
172
174 // Return EMCAL geometrical parameters
175 //
176
177 const Char_t* GetNameOfEMCALEnvelope() const { return "XEN1"; }
183 Float_t GetEnvelop(Int_t index) const { return mEnvelop[index]; }
185 Float_t GetZLength() const { return mZLength; }
191 Int_t GetNECLayers() const { return mNECLayers; }
193
196 Int_t GetNZ() const { return mNZ; }
197
200 Int_t GetNEta() const { return mNZ; }
201
204 Int_t GetNPhi() const { return mNPhi; }
205
208 Float_t GetSampling() const { return mSampling; }
219 Int_t GetNPhiSuperModule() const { return mNPhiSuperModule; }
220 Int_t GetNPHIdiv() const { return mNPHIdiv; }
221 Int_t GetNETAdiv() const { return mNETAdiv; }
222 Int_t GetNCells() const { return mNCells; }
224 Float_t GetTrd1Angle() const { return mTrd1Angle; }
225 Float_t Get2Trd1Dx2() const { return m2Trd1Dx2; }
228 // --
229 Int_t GetNCellsInSupMod() const { return mNCellsInSupMod; }
230 Int_t GetNCellsInModule() const { return mNCellsInModule; }
231 Int_t GetKey110DEG() const { return mKey110DEG; }
232 Int_t GetnSupModInDCAL() const { return mnSupModInDCAL; }
233 Int_t GetILOSS() const { return mILOSS; }
234 Int_t GetIHADR() const { return mIHADR; }
235 // --
238 Int_t GetNTowers() const { return mNPhi * mNZ; }
239 //
240 Double_t GetPhiCenterOfSM(Int_t nsupmod) const;
241 Double_t GetPhiCenterOfSMSec(Int_t nsupmod) const;
242 Float_t GetSuperModulesPar(Int_t ipar) const { return mParSM[ipar]; }
243 //
244 EMCALSMType GetSMType(Int_t nSupMod) const
245 {
246 if (nSupMod >= mNumberOfSuperModules) {
248 }
249 return mEMCSMSystem[nSupMod];
250 }
251
255 Bool_t IsDCALSM(Int_t nSupMod) const;
256
260 Bool_t IsDCALExtSM(Int_t nSupMod) const;
261
262 // Methods needed for SM in extension, where center of SM != center of the SM-section.
263 // Used in AliEMCALv0 to calculate position.
264 std::tuple<double, double> GetPhiBoundariesOfSM(Int_t nSupMod) const;
265 std::tuple<double, double> GetPhiBoundariesOfSMGap(Int_t nPhiSec) const;
266
267 // Obsolete?
269
271 // Geometry data member setters
272 //
273 void SetNZ(Int_t nz) { mNZ = nz; }
274 void SetNPhi(Int_t nphi) { mNPhi = nphi; }
275 //
276 void SetSampling(Float_t samp) { mSampling = samp; }
277
279 // Global geometry methods
280 //
281
288 void GetGlobal(const Double_t* loc, Double_t* glob, int ind) const;
289
296 void GetGlobal(const TVector3& vloc, TVector3& vglob, int ind) const;
297
304 void GetGlobal(Int_t absId, Double_t glob[3]) const;
305
311 void GetGlobal(Int_t absId, TVector3& vglob) const;
312
314 // May 31, 2006; ALICE numbering scheme:
315 // see ALICE-INT-2003-038: ALICE Coordinate System and Software Numbering Convention
316 // All indexes are stared from zero now.
317 //
318 // abs id <-> indexes; Shish-kebab case, only TRD1 now.
319 // EMCAL -> Super Module -> module -> tower(or cell) - logic tree of EMCAL
320 //
321 //** Usual name of variable - Dec 18,2006 **
322 // nSupMod - index of super module (SM)
323 // nModule - index of module in SM
324 // nIphi - phi index of tower(cell) in module
325 // nIeta - eta index of tower(cell) in module
326 //
327 // Inside SM
328 // iphim - phi index of module in SM
329 // ietam - eta index of module in SM
330 //
331 // iphi - phi index of tower(cell) in SM
332 // ieta - eta index of tower(cell) in SM
333 //
334 // for a given tower index absId returns eta and phi of gravity center of tower.
335
341 std::tuple<double, double> EtaPhiFromIndex(Int_t absId) const;
342
349 int GetAbsCellIdFromEtaPhi(Double_t eta, Double_t phi) const;
350
355 std::tuple<int, int> GlobalRowColFromIndex(int cellID) const;
356
361 int GlobalCol(int cellID) const;
362
367 int GlobalRow(int cellID) const;
368
374 int GetCellAbsIDFromGlobalRowCol(int row, int col) const;
375
381 std::tuple<int, int, int> GetPositionInSupermoduleFromGlobalRowCol(int row, int col) const;
382
388 std::tuple<int, int, int, int> GetCellIndexFromGlobalRowCol(int row, int col) const;
389
395 int SuperModuleNumberFromEtaPhi(Double_t eta, Double_t phi) const;
396
405 int GetAbsCellId(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const;
406
410 Bool_t CheckAbsCellId(Int_t absId) const;
411
416 std::tuple<int, int, int, int> GetCellIndex(Int_t absId) const;
417
422 std::tuple<int, int> GetModulePhiEtaIndexInSModule(int supermoduleID, int moduleID) const;
423
430 std::tuple<int, int> GetCellPhiEtaIndexInSModule(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const;
431
445 std::tuple<int, int> ShiftOnlineToOfflineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const;
446
458 std::tuple<int, int> ShiftOfflineToOnlineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const;
459
463 Int_t GetSuperModuleNumber(Int_t absId) const;
464 Int_t GetNumberOfModuleInPhiDirection(Int_t nSupMod) const
465 {
466 if (GetSMType(nSupMod) == EMCAL_HALF) {
467 return mNPhi / 2;
468 } else if (GetSMType(nSupMod) == EMCAL_THIRD) {
469 return mNPhi / 3;
470 } else if (GetSMType(nSupMod) == DCAL_EXT) {
471 return mNPhi / 3;
472 } else {
473 return mNPhi;
474 }
475 }
476
486 std::tuple<int, int, int> GetModuleIndexesFromCellIndexesInSModule(int supermoduleID, int phiInSupermodule, int etaInSupermodule) const;
487
493 Int_t GetAbsCellIdFromCellIndexes(Int_t nSupMod, Int_t iphi, Int_t ieta) const;
494
503 math_utils::Point3D<double> RelPosCellInSModule(Int_t absId, Double_t distEf) const;
504
510
516 std::tuple<int, int, int> getOnlineID(int towerID);
517
524 std::tuple<bool, int, int> areAbsIDsFromSameTCard(int absId1, int absId2) const;
525
532 std::tuple<int, int> getLinkAssignment(int ddlID) const { return std::make_tuple(mCRORCID[ddlID], mCRORCLink[ddlID]); };
533
534 std::vector<EMCALSMType> GetEMCSystem() const { return mEMCSMSystem; } // EMC System, SM type list
535 // Local Coordinates of SM
536 std::vector<Double_t> GetCentersOfCellsEtaDir() const
537 {
539 } // size fNEta*fNETAdiv (for TRD1 only) (eta or z in SM, in cm)
540 std::vector<Double_t> GetCentersOfCellsXDir() const
541 {
542 return mCentersOfCellsXDir;
543 } // size fNEta*fNETAdiv (for TRD1 only) ( x in SM, in cm)
544 std::vector<Double_t> GetCentersOfCellsPhiDir() const
545 {
547 } // size fNPhi*fNPHIdiv (for TRD1 only) (phi or y in SM, in cm)
548 //
549 std::vector<Double_t> GetEtaCentersOfCells() const
550 {
551 return mEtaCentersOfCells;
552 } // [fNEta*fNETAdiv*fNPhi*fNPHIdiv], positive direction (eta>0); eta depend from phi position;
553 std::vector<Double_t> GetPhiCentersOfCells() const
554 {
555 return mPhiCentersOfCells;
556 } // [fNPhi*fNPHIdiv] from center of SM (-10. < phi < +10.)
557
559 // useful utilities
560 //
562 { // returns theta in radians for a given pseudorapidity
563 return 2.0 * TMath::ATan(TMath::Exp(-eta));
564 }
566 { // returns z in for a given
567 // pseudorapidity and r=sqrt(x*x+y*y).
568 return r / TMath::Tan(AngleFromEta(eta));
569 }
570
575 void SetMisalMatrix(const TGeoHMatrix* m, Int_t smod) const;
576
580 void SetMisalMatrixFromCcdb(const char* path = "Users/m/mhemmer/EMCAL/Config/GeometryAligned", int timestamp = 10000) const;
581
591 void RecalculateTowerPosition(Float_t drow, Float_t dcol, const Int_t sm, const Float_t depth,
592 const Float_t misaligTransShifts[15], const Float_t misaligRotShifts[15],
593 Float_t global[3]) const;
594
599 const TGeoHMatrix* GetMatrixForSuperModule(Int_t smod) const;
600
604 const TGeoHMatrix* GetMatrixForSuperModuleFromGeoManager(Int_t smod) const;
605
612 const TGeoHMatrix* GetMatrixForSuperModuleFromArray(Int_t smod) const;
613
614 protected:
616 void Init();
617
619 void DefineEMC(std::string_view mcname, std::string_view mctitle);
620
627 std::tuple<int, int, int, int> CalculateCellIndex(Int_t absId) const;
628
629 std::string mGeoName;
633 Int_t mNETAdiv;
634 Int_t mNPHIdiv;
636 std::vector<Double_t> mPhiBoundariesOfSM;
637 std::vector<Double_t> mPhiCentersOfSM;
638 std::vector<Double_t> mPhiCentersOfSMSec;
639
640 // Local Coordinates of SM
641 std::vector<Double_t> mPhiCentersOfCells;
642 std::vector<Double_t> mCentersOfCellsEtaDir;
643 std::vector<Double_t> mCentersOfCellsPhiDir;
644 std::vector<Double_t>
646 Int_t mNCells;
647 Int_t mNPhi;
648 std::vector<Double_t> mCentersOfCellsXDir;
661 std::vector<ShishKebabTrd1Module> mShishKebabTrd1Modules;
667 Int_t mNZ;
670
671 // Geometry Parameters
675
676 // Members from the EMCGeometry class
680
681 // Shish-kebab option - 23-aug-04 by PAI; COMPACT, TWIST, TRD1 and TRD2
683
685 std::vector<EMCALSMType> mEMCSMSystem;
686
690
693
694 // TRD1 options - 30-sep-04
698
699 // Oct 26,2010
702
703 Int_t mILOSS;
704 Int_t mIHADR;
705
707
708 std::array<int, 46> mCRORCID = {110, 110, 112, 112, 110, 110, 112, 112, 110, 110, 112, 112, 111, 111, 113, 113, 111, 111, 113, 113, 111, 111, 113, 113, 114, 114, 116, 116, 114, 114, 116, 116, 115, 115, 117, 117, 115, 115, 117, 117, -1, -1, -1, -1, 111, 117}; // CRORC ID w.r.t SM
709 std::array<int, 46> mCRORCLink = {0, 1, 0, 1, 2, 3, 2, 3, 4, 5, 4, 5, 0, 1, 0, 1, 2, 3, 2, 3, 4, -1, 4, 5, 0, 1, 0, 1, 2, 3, 2, 3, 0, 1, 0, 1, 2, 3, 2, -1, -1, -1, -1, -1, 5, 3}; // CRORC limk w.r.t FEE ID
710
711 mutable const TGeoHMatrix* SMODULEMATRIX[EMCAL_MODULES];
712 std::vector<std::tuple<int, int, int, int>> mCellIndexLookup;
713
714 private:
715 static Geometry* sGeom;
716};
717
718inline Bool_t Geometry::CheckAbsCellId(Int_t absId) const
719{
720 if (absId < 0 || absId >= mNCells) {
721 return kFALSE;
722 } else {
723 return kTRUE;
724 }
725}
726} // namespace emcal
727} // namespace o2
728#endif
uint32_t col
Definition RawData.h:4
EMCAL geometry definition.
Definition Geometry.h:42
Float_t GetDCALInnerExtandedEta() const
Definition Geometry.h:192
std::tuple< int, int, int > GetModuleIndexesFromCellIndexesInSModule(int supermoduleID, int phiInSupermodule, int etaInSupermodule) const
Transition from cell indexes (iphi, ieta) to module indexes (iphim, ietam, nModule)
Definition Geometry.cxx:801
Float_t mFrontSteelStrip
13-may-05
Definition Geometry.h:687
static Geometry * GetInstanceFromRunNumber(Int_t runNumber, const std::string_view="", const std::string_view mcname="TGeant3", const std::string_view mctitle="")
Instanciate geometry depending on the run number. Mostly used in analysis and MC anchors.
Definition Geometry.cxx:242
Float_t mEtaMaxOfTRD1
Max eta in case of TRD1 geometry (see AliEMCALShishKebabTrd1Module)
Definition Geometry.h:654
Int_t GetSuperModuleNumber(Int_t absId) const
Get cell SM, from absolute ID number.
void RecalculateTowerPosition(Float_t drow, Float_t dcol, const Int_t sm, const Float_t depth, const Float_t misaligTransShifts[15], const Float_t misaligRotShifts[15], Float_t global[3]) const
Int_t GetNCells() const
Definition Geometry.h:222
std::tuple< int, int, int, int > CalculateCellIndex(Int_t absId) const
Calculate cell SM, module numbers from absolute ID number.
Float_t mIPDistance
Radial Distance of the inner surface of the EMCAL.
Definition Geometry.h:668
Float_t mDCALPhiMin
Minimum angular position of DCAL in Phi (degrees)
Definition Geometry.h:655
Float_t mPassiveScintThick
13-may-05
Definition Geometry.h:689
Float_t mPhiTileSize
Size of phi tile.
Definition Geometry.h:665
int GlobalRow(int cellID) const
Get row number of cell in global numbering scheme.
Definition Geometry.cxx:919
Float_t GetDCALPhiMax() const
Definition Geometry.h:188
Int_t mNZ
Number of Towers in the Z direction.
Definition Geometry.h:667
Float_t mECScintThick
cm, Thickness of the scintillators
Definition Geometry.h:678
const TGeoHMatrix * GetMatrixForSuperModuleFromGeoManager(Int_t smod) const
Provides shift-rotation matrix for EMCAL from the TGeoManager.
std::tuple< int, int > GetCellPhiEtaIndexInSModule(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const
Get eta-phi indexes of cell in SM.
Float_t GetSampling() const
Definition Geometry.h:208
std::tuple< int, int, int, int > GetCellIndex(Int_t absId) const
Get cell SM, module numbers from absolute ID number.
Float_t GetArm1PhiMin() const
Definition Geometry.h:178
Int_t GetNCellsInSupMod() const
Definition Geometry.h:229
std::tuple< double, double > GetPhiBoundariesOfSM(Int_t nSupMod) const
Int_t mNECLayers
number of scintillator layers
Definition Geometry.h:679
Int_t GetNZ() const
Get the number of modules in supermodule in z- (beam) direction.
Definition Geometry.h:196
std::vector< Double_t > mPhiCentersOfCells
[fNPhi*fNPHIdiv] from center of SM (-10. < phi < +10.)
Definition Geometry.h:641
const TGeoHMatrix * GetMatrixForSuperModule(Int_t smod) const
Provides shift-rotation matrix for EMCAL from externally set matrix or from TGeoManager.
Geometry & operator=(const Geometry &rvalue)
Assignment operator.
Definition Geometry.cxx:193
Float_t GetDCALInnerEdge() const
Definition Geometry.h:186
static Bool_t IsInitialized()
Definition Geometry.h:115
std::tuple< int, int, int > getOnlineID(int towerID)
Get link ID, row and column from cell ID, have a look here: https://alice.its.cern....
Int_t GetNEta() const
Get the number of modules in supermodule in #eta direction.
Definition Geometry.h:200
Float_t GetLateralSteelStrip() const
Definition Geometry.h:214
void SetSampling(Float_t samp)
Definition Geometry.h:276
Float_t mEtaModuleSize
Eta -> Y.
Definition Geometry.h:664
void SetMisalMatrixFromCcdb(const char *path="Users/m/mhemmer/EMCAL/Config/GeometryAligned", int timestamp=10000) const
Double_t GetPhiCenterOfSMSec(Int_t nsupmod) const
Float_t GetPhiModuleSize() const
Definition Geometry.h:211
Float_t GetSuperModulesPar(Int_t ipar) const
Definition Geometry.h:242
std::vector< Double_t > GetCentersOfCellsXDir() const
Definition Geometry.h:540
Bool_t CheckAbsCellId(Int_t absId) const
Check whether a cell number is valid.
Definition Geometry.h:718
Float_t Get2Trd1Dx2() const
Definition Geometry.h:225
const Char_t * GetNameOfEMCALEnvelope() const
Definition Geometry.h:177
Float_t GetFrontSteelStrip() const
Definition Geometry.h:213
Float_t mZLength
Total length in z direction.
Definition Geometry.h:673
Float_t mEnvelop[3]
The GEANT TUB for the detector.
Definition Geometry.h:649
Float_t mDCALInnerEdge
Inner edge for DCAL.
Definition Geometry.h:660
void CreateListOfTrd1Modules()
std::vector< Double_t > GetCentersOfCellsPhiDir() const
Definition Geometry.h:544
Float_t mParSM[3]
SM sizes as in GEANT (TRD1)
Definition Geometry.h:662
Float_t GetShellThickness() const
Definition Geometry.h:184
Int_t mNCellsInModule
Number cell in module.
Definition Geometry.h:635
std::vector< EMCALSMType > mEMCSMSystem
geometry structure
Definition Geometry.h:685
Bool_t IsDCALExtSM(Int_t nSupMod) const
Check if iSupMod is a valid DCal 1/3rd SM.
std::tuple< int, int, int, int > GetCellIndexFromGlobalRowCol(int row, int col) const
Get the cell indices from global position in the EMCAL.
Definition Geometry.cxx:903
Int_t mIHADR
Options for Geant (MIP business) - will call in AliEMCAL.
Definition Geometry.h:704
Double_t GetPhiCenterOfSM(Int_t nsupmod) const
Int_t GetKey110DEG() const
Definition Geometry.h:231
const TGeoHMatrix * SMODULEMATRIX[EMCAL_MODULES]
Orientations of EMCAL super modules.
Definition Geometry.h:711
std::tuple< int, int > GetModulePhiEtaIndexInSModule(int supermoduleID, int moduleID) const
Get eta-phi indexes of module in SM.
std::tuple< int, int, int > GetPositionInSupermoduleFromGlobalRowCol(int row, int col) const
Get the posision (row, col) of a global row-col position.
Definition Geometry.cxx:857
Float_t GetArm1EtaMin() const
Definition Geometry.h:180
const TGeoHMatrix * GetMatrixForSuperModuleFromArray(Int_t smod) const
Provides shift-rotation matrix for EMCAL from fkSModuleMatrix[smod].
Float_t GetPhiSuperModule() const
Definition Geometry.h:218
Float_t ZFromEtaR(Float_t r, Float_t eta) const
Definition Geometry.h:565
Float_t mTrd1AlFrontThick
Thickness of the Al front plate.
Definition Geometry.h:700
Float_t mSampling
Sampling factor.
Definition Geometry.h:674
Float_t mDCALPhiMax
Maximum angular position of DCAL in Phi (degrees)
Definition Geometry.h:656
Float_t mTrd1BondPaperThick
Thickness of the Bond Paper sheet.
Definition Geometry.h:701
Int_t GetNCellsInModule() const
Definition Geometry.h:230
Float_t GetPhiGapForSuperModules() const
Definition Geometry.h:210
Float_t GetEMCALPhiMax() const
Definition Geometry.h:189
Int_t mNCells
Number of cells in calo.
Definition Geometry.h:646
std::tuple< bool, int, int > areAbsIDsFromSameTCard(int absId1, int absId2) const
Check if 2 cells belong to the same T-Card.
Float_t mArm1PhiMin
Minimum angular position of EMCAL in Phi (degrees)
Definition Geometry.h:652
void SetMisalMatrix(const TGeoHMatrix *m, Int_t smod) const
Float_t GetDeltaEta() const
Definition Geometry.h:236
const std::vector< ShishKebabTrd1Module > & GetShishKebabTrd1Modules() const
Definition Geometry.h:124
Int_t mNETAdiv
Number eta division of module.
Definition Geometry.h:633
std::vector< EMCALSMType > GetEMCSystem() const
Definition Geometry.h:534
Float_t mLateralSteelStrip
13-may-05
Definition Geometry.h:688
std::tuple< int, int > GlobalRowColFromIndex(int cellID) const
get (Column,Row) pair of cell in global numbering scheme
Definition Geometry.cxx:831
Float_t GetDCALPhiMin() const
Definition Geometry.h:187
Int_t GetnSupModInDCAL() const
Definition Geometry.h:232
Float_t mTrd1Angle
angle in x-z plane (in degree)
Definition Geometry.h:695
Float_t mECPbRadThickness
cm, Thickness of the Pb radiators
Definition Geometry.h:677
void GetGlobal(const Double_t *loc, Double_t *glob, int ind) const
Figure out the global coordinates from local coordinates on a supermodule.
Definition Geometry.cxx:704
std::vector< Double_t > mEtaCentersOfCells
[fNEta*fNETAdiv*fNPhi*fNPHIdiv], positive direction (eta>0); eta depend from phi position;
Definition Geometry.h:645
int GetAbsCellId(int supermoduleID, int moduleID, int phiInModule, int etaInModule) const
Get cell absolute ID number from location module (2 times 2 cells) of a super module.
Definition Geometry.cxx:761
Float_t mArm1PhiMax
Maximum angular position of EMCAL in Phi (degrees)
Definition Geometry.h:653
static Geometry * GetInstance()
Get geometry instance. It should have been set before.
Definition Geometry.cxx:213
std::array< int, 46 > mCRORCID
Definition Geometry.h:708
std::tuple< double, double > EtaPhiFromIndex(Int_t absId) const
Figure out the eta/phi coordinates of a cell.
Definition Geometry.cxx:754
Int_t mKey110DEG
For calculation abs cell id; 19-oct-05.
Definition Geometry.h:630
std::vector< Double_t > mPhiCentersOfSMSec
Phi of centers of section where SM lies; size is fNumberOfSuperModules/2.
Definition Geometry.h:638
std::vector< Double_t > mCentersOfCellsEtaDir
Size fNEta*fNETAdiv (for TRD1 only) (eta or z in SM, in cm)
Definition Geometry.h:642
std::string mGeoName
Geometry name string.
Definition Geometry.h:629
Int_t GetNPhiSuperModule() const
Definition Geometry.h:219
int GlobalCol(int cellID) const
Get column number of cell in global numbering scheme.
Definition Geometry.cxx:914
Float_t mShellThickness
Total thickness in (x,y) direction.
Definition Geometry.h:672
EMCALSMType GetSMType(Int_t nSupMod) const
Definition Geometry.h:244
std::vector< Double_t > mCentersOfCellsPhiDir
Size fNPhi*fNPHIdiv (for TRD1 only) (phi or y in SM, in cm)
Definition Geometry.h:643
std::vector< Double_t > GetPhiCentersOfCells() const
Definition Geometry.h:553
Int_t GetILOSS() const
Definition Geometry.h:233
Float_t GetSteelFrontThickness() const
Definition Geometry.h:268
Int_t mNPHIdiv
Number phi division of module.
Definition Geometry.h:634
Int_t mnSupModInDCAL
For calculation abs cell id; 06-nov-12.
Definition Geometry.h:631
Float_t mEtaTileSize
Size of eta tile.
Definition Geometry.h:666
std::tuple< int, int > getLinkAssignment(int ddlID) const
Temporary link assignment (till final link assignment is known -.
Definition Geometry.h:532
Float_t GetDeltaPhi() const
Definition Geometry.h:237
Float_t GetLongModuleSize() const
Definition Geometry.h:223
std::array< int, 46 > mCRORCLink
Definition Geometry.h:709
Float_t GetPhiTileSize() const
Definition Geometry.h:216
int SuperModuleNumberFromEtaPhi(Double_t eta, Double_t phi) const
Given a global eta/phi point check if it belongs to a supermodule covered region.
Definition Geometry.cxx:924
void DefineEMC(std::string_view mcname, std::string_view mctitle)
Init function of previous class EMCGeometry.
Definition Geometry.cxx:366
Bool_t IsDCALSM(Int_t nSupMod) const
Check if iSupMod is a valid DCal standard SM.
void Init()
initializes the parameters of EMCAL
Int_t GetNTowers() const
Definition Geometry.h:238
Bool_t Impact(const TParticle *particle) const
Check if particle falls in the EMCal/DCal geometry.
std::vector< ShishKebabTrd1Module > mShishKebabTrd1Modules
List of modules.
Definition Geometry.h:661
Float_t GetZLength() const
Definition Geometry.h:185
Bool_t IsInEMCAL(const math_utils::Point3D< double > &pnt) const
Checks whether point is inside the EMCal volume.
Float_t AngleFromEta(Float_t eta) const
Definition Geometry.h:561
math_utils::Point3D< double > RelPosCellInSModule(Int_t absId, Double_t distEf) const
Look to see what the relative position inside a given cell is for a recpoint.
void SetNZ(Int_t nz)
Definition Geometry.h:273
const std::string & GetName() const
Definition Geometry.h:111
Int_t GetIHADR() const
Definition Geometry.h:234
Float_t mEMCALPhiMax
Maximum angular position of EMCAL in Phi (degrees)
Definition Geometry.h:657
Bool_t IsInDCAL(const math_utils::Point3D< double > &pnt) const
Checks whether point is inside the DCal volume.
static const std::string & GetDefaultGeometryName()
Definition Geometry.h:113
Int_t GetNECLayers() const
Definition Geometry.h:191
Float_t GetIPDistance() const
Definition Geometry.h:182
Int_t mNPhiSuperModule
9 - number supermodule in phi direction
Definition Geometry.h:692
Float_t mPhiModuleSize
Phi -> X.
Definition Geometry.h:663
Float_t mPhiGapForSM
Gap betweeen supermodules in phi direction.
Definition Geometry.h:697
Geometry()=default
Default constructor. It must be kept public for root persistency purposes, but should never be called...
Float_t m2Trd1Dx2
2*dx2 for TRD1
Definition Geometry.h:696
Float_t GetDCALStandardPhiMax() const
Definition Geometry.h:190
AcceptanceType_t IsInEMCALOrDCAL(const math_utils::Point3D< double > &pnt) const
Checks whether point is inside the EMCal volume (included DCal)
Float_t GetTrd1Angle() const
Definition Geometry.h:224
std::vector< Double_t > GetEtaCentersOfCells() const
Definition Geometry.h:549
Int_t mILOSS
Options for Geant (MIP business) - will call in AliEMCAL.
Definition Geometry.h:703
Int_t GetNumberOfSuperModules() const
Definition Geometry.h:209
Float_t GetECScintThick() const
Definition Geometry.h:207
~Geometry()
Destructor.
Definition Geometry.cxx:199
Int_t GetNumberOfModuleInPhiDirection(Int_t nSupMod) const
Definition Geometry.h:464
std::tuple< int, int > ShiftOnlineToOfflineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const
Adapt cell indices in supermodule to online indexing.
Float_t mPhiSuperModule
Phi of normal supermodule (20, in degree)
Definition Geometry.h:691
Float_t GetTrd1BondPaperThick() const
Definition Geometry.h:227
Float_t mArm1EtaMin
Minimum pseudorapidity position of EMCAL in Eta.
Definition Geometry.h:650
void ImpactOnEmcal(const math_utils::Point3D< double > &vtx, Double_t theta, Double_t phi, Int_t &absId, math_utils::Point3D< double > &vimpact) const
Get the impact coordinates on EMCAL.
Float_t GetArm1EtaMax() const
Definition Geometry.h:181
std::vector< Double_t > GetCentersOfCellsEtaDir() const
Definition Geometry.h:536
Int_t GetNPhi() const
Get the number of modules in supermodule in #phi direction.
Definition Geometry.h:204
Float_t mArm1EtaMax
Maximum pseudorapidity position of EMCAL in Eta.
Definition Geometry.h:651
Int_t GetAbsCellIdFromCellIndexes(Int_t nSupMod, Int_t iphi, Int_t ieta) const
Transition from super module number (nSupMod) and cell indexes (ieta,iphi) to cell absolute ID number...
Definition Geometry.cxx:814
std::vector< Double_t > mCentersOfCellsXDir
Size fNEta*fNETAdiv (for TRD1 only) ( x in SM, in cm)
Definition Geometry.h:648
void SetNPhi(Int_t nphi)
Definition Geometry.h:274
Float_t GetEtaModuleSize() const
Definition Geometry.h:212
Float_t GetTrd1AlFrontThick() const
Definition Geometry.h:226
Int_t mNPhi
Number of Towers in the PHI direction.
Definition Geometry.h:647
const ShishKebabTrd1Module & GetShishKebabModule(Int_t neta) const
Get the Module parameters for a eta.
Int_t mNumberOfSuperModules
default is 12 = 6 * 2
Definition Geometry.h:682
std::tuple< double, double > GetPhiBoundariesOfSMGap(Int_t nPhiSec) const
std::vector< Double_t > mPhiBoundariesOfSM
Phi boundaries of SM in rad; size is fNumberOfSuperModules;.
Definition Geometry.h:636
std::vector< std::tuple< int, int, int, int > > mCellIndexLookup
Lookup table for cell indices.
Definition Geometry.h:712
void DefineSamplingFraction(const std::string_view mcname="", const std::string_view mctitle="")
Set the value of the Sampling used to calibrate the MC hits energy (check)
Definition Geometry.cxx:308
Int_t GetNPHIdiv() const
Definition Geometry.h:220
Float_t mDCALStandardPhiMax
Special edge for the case that DCAL contian extension.
Definition Geometry.h:658
Float_t mLongModuleSize
Size of long module.
Definition Geometry.h:669
Float_t GetEnvelop(Int_t index) const
Definition Geometry.h:183
Float_t mSteelFrontThick
Thickness of the front stell face of the support box - 9-sep-04; obsolete?
Definition Geometry.h:706
Float_t GetECPbRadThick() const
Definition Geometry.h:206
Float_t GetArm1PhiMax() const
Definition Geometry.h:179
int GetAbsCellIdFromEtaPhi(Double_t eta, Double_t phi) const
Get cell absolute ID number from eta and phi location.
Definition Geometry.cxx:955
Int_t mNCellsInSupMod
Number cell in super module.
Definition Geometry.h:632
int GetCellAbsIDFromGlobalRowCol(int row, int col) const
Get the absolute cell ID from global position in the EMCAL.
Definition Geometry.cxx:897
std::tuple< int, int > ShiftOfflineToOnlineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const
Adapt cell indices in supermodule to offline indexing.
std::vector< Double_t > mPhiCentersOfSM
Phi of centers of SM; size is fNumberOfSuperModules/2.
Definition Geometry.h:637
Float_t GetPassiveScintThick() const
Definition Geometry.h:215
Float_t mDCALInnerExtandedEta
DCAL inner edge in Eta (with some extension)
Definition Geometry.h:659
Float_t GetEtaTileSize() const
Definition Geometry.h:217
Int_t GetNETAdiv() const
Definition Geometry.h:221
Main class for TRD1 geometry of Shish-Kebab case.
Handling error due to invalid supermodule.
const GLfloat * m
Definition glcorearb.h:4066
GLuint index
Definition glcorearb.h:781
GLuint const GLchar * name
Definition glcorearb.h:781
GLint GLint GLsizei GLsizei GLsizei depth
Definition glcorearb.h:470
GLsizei const GLchar *const * path
Definition glcorearb.h:3591
GLboolean r
Definition glcorearb.h:1233
const std::string DEFAULT_GEOMETRY
@ EMCAL_MODULES
Number of modules, 12 for EMCal + 8 for DCAL.
Definition Constants.h:24
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::vector< int > row