ATLAS Offline Software
Loading...
Searching...
No Matches
HVHelper.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef LARG4_EC_HVHELPER_H
6#define LARG4_EC_HVHELPER_H
7
8#include <stdio.h>
9#include <stdexcept>
10
12#include "AthenaKernel/Units.h"
13
15
16#include "G4ThreeVector.hh"
17#include "globals.hh"
18
19#include <memory>
20
21namespace LArG4 { namespace EC {
22
23class HVHelper : public AthMessaging
24{
25 protected:
26 static constexpr G4int s_NofAtlasSide = 2;
27 static constexpr G4int s_NofEtaSection = 9;
28 static constexpr G4int s_NofElectrodeSide = 2;
29 static constexpr G4int s_NofElectrodesOut = 768;
30
31 static const G4double s_EtaLimit[s_NofEtaSection + 1];
32
34
36 const LArWheelCalculator *lwc(void) const { return m_calculator; }
37
38 const G4double m_WheelShift;
39 const G4int m_NofPhiSections;
40
41 const G4int m_EtaMin, m_EtaMax;
42
44 const LArWheelCalculator *calc,
45 const G4String &version
46 );
47 void AcquireMaps(const G4String &version, G4bool from_DB);
48 void GetMapFromDB(void);
49 virtual void ReadMapFromFile(const G4String &version) = 0;
50 FILE * OpenFileAndCheckVersion(const G4String &version);
51
52 G4int GetEtaSection(const G4ThreeVector &p) const;
53
54 virtual G4int GetPhiSection(
55 G4int, G4int, G4int, G4int
56 ) const = 0;
57
58 public:
60 G4double phi, G4int compartment, G4int eta_bin
61 ) const;
62
63 G4double GetVoltage(const G4ThreeVector& p) const;
64
65 virtual G4double GetVoltage(const G4ThreeVector&, G4int, G4int) const;
66 virtual G4double GetVoltage(
67 const G4ThreeVector&, const std::pair<G4int, G4int> &
68 ) const;
69
70 static std::unique_ptr<const HVHelper> CreateHelper(
71 const LArWheelCalculator *calc,
72 const G4String &version,
73 G4bool fromDB
74 );
75
76 virtual ~HVHelper() {}
77};
78
79class HVHelperV00 : public HVHelper // This serves map versions 00 and 01
80{
81 private:
83 G4int StartPhi(G4int side, G4int eta, G4int ele) const
84 {
85 return m_StartPhi[(side*s_NofEtaSection + eta)*s_NofElectrodeSide + ele];
86 }
87
88 void SetStartPhi(G4int value, G4int side, G4int eta, G4int ele)
89 {
90 m_StartPhi[(side*s_NofEtaSection + eta)*s_NofElectrodeSide + ele] = value;
91 }
92
94
95 void ReadMapFromFile(const G4String &version) override final;
96
97 virtual G4int GetPhiSection(
98 G4int, G4int, G4int, G4int
99 ) const override final;
100
101 public:
103 const LArWheelCalculator *calc,
104 const G4String &version,
105 G4bool fromDB
106 ) : HVHelper(calc, version)
108 {
109 AcquireMaps(version, fromDB);
110 }
111};
112
113class HVHelperV02 : public HVHelper // This serves map versions 02 and above
114{
115 private:
116 void ReadMapFromFile(const G4String &version) override final;
117
118 virtual G4int GetPhiSection(
119 G4int, G4int, G4int, G4int
120 ) const override final;
121
122 public:
124 const LArWheelCalculator *calc,
125 const G4String &version,
126 G4bool fromDB
127 ) : HVHelper(calc, version) { AcquireMaps(version, fromDB); }
128};
129
130}} // namespaces
131
132#endif // LARG4_EC_HVHELPER_H
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Wrapper to avoid constant divisions when using units.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
G4int m_StartPhi[s_NofAtlasSide *s_NofEtaSection *s_NofElectrodeSide]
Definition HVHelper.h:82
G4int StartPhi(G4int side, G4int eta, G4int ele) const
Definition HVHelper.h:83
const G4int m_NumberOfElectrodesInPhiSection
Definition HVHelper.h:93
void SetStartPhi(G4int value, G4int side, G4int eta, G4int ele)
Definition HVHelper.h:88
HVHelperV00(const LArWheelCalculator *calc, const G4String &version, G4bool fromDB)
Definition HVHelper.h:102
virtual G4int GetPhiSection(G4int, G4int, G4int, G4int) const override final
void ReadMapFromFile(const G4String &version) override final
HVHelperV02(const LArWheelCalculator *calc, const G4String &version, G4bool fromDB)
Definition HVHelper.h:123
virtual G4int GetPhiSection(G4int, G4int, G4int, G4int) const override final
void ReadMapFromFile(const G4String &version) override final
static constexpr G4int s_NofElectrodesOut
Definition HVHelper.h:29
const G4int m_NofPhiSections
Definition HVHelper.h:39
G4double m_Values[s_NofAtlasSide][s_NofEtaSection][s_NofElectrodeSide][s_NofElectrodesOut]
Definition HVHelper.h:33
HVHelper(const LArWheelCalculator *calc, const G4String &version)
static constexpr G4int s_NofEtaSection
Definition HVHelper.h:27
static constexpr G4int s_NofAtlasSide
Definition HVHelper.h:26
const LArWheelCalculator * lwc(void) const
Definition HVHelper.h:36
virtual G4double GetVoltage(const G4ThreeVector &, G4int, G4int) const
static std::unique_ptr< const HVHelper > CreateHelper(const LArWheelCalculator *calc, const G4String &version, G4bool fromDB)
static constexpr G4int s_NofElectrodeSide
Definition HVHelper.h:28
virtual ~HVHelper()
Definition HVHelper.h:76
virtual void ReadMapFromFile(const G4String &version)=0
void GetMapFromDB(void)
G4double GetVoltage(const G4ThreeVector &p) const
const G4int m_EtaMin
Definition HVHelper.h:41
static const G4double s_EtaLimit[s_NofEtaSection+1]
Definition HVHelper.h:31
virtual G4double GetVoltage(const G4ThreeVector &, const std::pair< G4int, G4int > &) const
FILE * OpenFileAndCheckVersion(const G4String &version)
const G4double m_WheelShift
Definition HVHelper.h:38
G4double GetVoltageBarrett(G4double phi, G4int compartment, G4int eta_bin) const
virtual G4int GetPhiSection(G4int, G4int, G4int, G4int) const =0
const G4int m_EtaMax
Definition HVHelper.h:41
G4int GetEtaSection(const G4ThreeVector &p) const
void AcquireMaps(const G4String &version, G4bool from_DB)
const LArWheelCalculator * m_calculator
Definition HVHelper.h:35
This class separates some of the geometry details of the LAr endcap.