ATLAS Offline Software
Loading...
Searching...
No Matches
EnergyCalculator.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// EnergyCalculator.h
6// Prepared 10-Apr-2002 Bill Seligman
7// from code written by Jozsef Toth.
8// 07-May-2003 AMS: now EnergyCalculator is not a singleton
9// 02-July-2003 J.T. Charge collection added
10// 20-July-2004 J.T. FieldMapVersion variable is added
11// 25-May-2005 J.T. -calling sequence of GapAdj changed
12// -new variables and function for reading and handling
13// values of HV of power supplies
14// - IonReco :for suppress signal because of ion recombination
15// - DriftVelo: Walkowiak's formula for drift velocity
16// Sept-2006 J.T. - collect signal from the Barrette volume
17// Nov -2006 J.T. - fieldmap array structures changed,
18// - lengths defined dinamically
19// - field map for first/last fold and for Barrett volume are included
20// May 2009 AMS move to namespace LArG4::EC
21// duplicated data members removed
22
23#ifndef LArG4_EC_EnergyCalculator_H
24#define LArG4_EC_EnergyCalculator_H
25
26#include <string>
27#include <stdexcept>
28
29#include "CLHEP/Units/SystemOfUnits.h"
30
35
38
41
42#include "G4ThreeVector.hh"
43#include "HVHelper.h"
44#include "globals.hh"
45
47class LArG4BirksLaw;
48
49namespace LArG4 {
50
51 namespace EC {
52
54 {
55 public:
56
57 EnergyCalculator(const std::string& name, ISvcLocator *pSvcLocator);
58 // Update handlers
59 void CorrectionTypeHandler(Gaudi::Details::PropertyBase&);
60 void SolidTypeHandler(Gaudi::Details::PropertyBase&);
61
62 virtual StatusCode initialize() override final;
63 virtual StatusCode finalize() override final;
64
66 // The interface for LArVCalculator.
67 virtual G4float OOTcut() const override final { return m_OOTcut; }
68 virtual G4bool Process(const G4Step*, std::vector<LArHitData>&) const override final;
69
70 // Check if the current hitTime is in-time
71 virtual G4bool isInTime(G4double hitTime) const override final
72 {
73 return !(hitTime > m_OOTcut); //FIXME should we be checking the absolute value of hitTime here?
74 }
75
76
77 private:
78 G4bool (EnergyCalculator::*m_Process_type) (const G4Step*, std::vector<LArHitData>&) const{};
79 G4double (EnergyCalculator::*m_GetGapSize_type) (const G4ThreeVector &p) const{};
80 G4double (EnergyCalculator::*m_distance_to_the_nearest_electrode_type) (const G4ThreeVector &p, G4double /*Barret_PhiStart*/) const{};
81
82 G4bool Process_Default(const G4Step*, std::vector<LArHitData>&) const;
83 G4bool Process_Barrett(const G4Step*, std::vector<LArHitData>&) const;
84 G4bool FindIdentifier_Default(const G4Step *, std::vector<LArHitData>&, G4ThreeVector &, G4ThreeVector &) const;
85 G4bool FindIdentifier_Barrett(const G4Step *, G4double, std::vector<LArHitData>&, G4ThreeVector &, G4ThreeVector &) const;
86 G4bool FindDMIdentifier_Barrett(const G4Step* step, std::vector<LArHitData>&) const;
87 G4bool GetCompartment_Barrett(const G4ThreeVector&,G4double,G4double,G4double,
88 G4int &, G4int &) const;
89
90 G4double GetHV_Value(const G4ThreeVector& p, G4double PhiStartOfPhiDiv) const;
91 G4double GetGapSize_Default(const G4ThreeVector &p) const {
92 return GetGapSize(p);
93 }
94 G4double GetGapSize_Barrett(const G4ThreeVector &p) const;
95 G4int GetPhiGap_Barrett(const G4ThreeVector &p, G4double PhiStartOfPhiDiv) const;
96 G4double distance_to_the_nearest_electrode_Default(const G4ThreeVector &p, G4double /*Barret_PhiStart*/) const {
98 }
99 G4double distance_to_the_nearest_electrode_Barrett(const G4ThreeVector &p, G4double Barret_PhiStart) const;
100
101 ServiceHandle<ILArCalibCalculatorSvc> m_supportCalculator{this, "SupportCalculator", "EMECSupportCalibrationCalculator"};
102
103 void SetConst_OuterBarrett(void); // used only for initialization
104 G4bool GetVolumeIndex(const G4Step *, G4int &, G4int &) const;
105 static const G4double s_LongBarThickness;// = 20. *mm;
106 static const G4double s_ColdCorrection;// =1.0036256;
107 static const G4double s_StripWidth;// =3.*mm/ColdCorrection;
108 static const G4double s_KapGap;// =1.*mm/ColdCorrection;
109 static const G4double s_EdgeWidth;// =1.*mm;
110 static const G4double s_DistOfEndofCuFromBack;// =22.77*mm/ColdCorrection;
111 static const G4double s_DistOfStartofCuFromBack;//=31.*mm; // frontface of the barrette
112 static const G4double s_ZmaxOfSignal;// DistOfStartofCuFromBack - DistOfEndofCuFromBack + EdgeWidth;
113
114 static const G4double s_S3_Etalim[21];
115 static const G4double s_Rmeas_outer[50];
116 static const G4double s_Zmeas_outer[2];
117 G4double m_RefzDist = 0.0; // = dElecFocaltoWRP+dWRPtoFrontFace+WheelThickness+ // used as const after
118 // +dWRPtoFrontFace+ LongBarThickness // initialization
119 // -DistOfEndofCuFromBack
120
121 G4double m_S3_Rlim[21]{}; // used as const after init
122 G4double m_rlim[50]{}; // used as const after init
123 G4double m_zlim[4]{}; // used as const after init
124
125 UnsignedIntegerProperty m_corrProp{this, "EnergyCorrection", 8, &EnergyCalculator::CorrectionTypeHandler};
127
128 G4double (EnergyCalculator::*m_ecorr_method) (G4double, const G4ThreeVector&, const G4ThreeVector&, G4double /*Barret_PhiStart*/) const{};
129 G4double dummy_correction_method(G4double e, const G4ThreeVector&, const G4ThreeVector&,
130 G4double /*Barret_PhiStart*/) const {
131 return e;
132 }
133 G4double GapAdjustment_old(G4double, const G4ThreeVector&, const G4ThreeVector&, G4double /*Barret_PhiStart*/) const;
134 G4double GapAdjustment (G4double, const G4ThreeVector&, const G4ThreeVector&, G4double /*Barret_PhiStart*/) const;
135 G4double GapAdjustment_E (G4double, const G4ThreeVector&, const G4ThreeVector&, G4double /*Barret_PhiStart*/) const;
136 G4double GapAdjustment_s (G4double, const G4ThreeVector&, const G4ThreeVector&, G4double /*Barret_PhiStart*/) const;
137 G4double GapAdjustment__sE (G4double, const G4ThreeVector&, const G4ThreeVector&, G4double /*Barret_PhiStart*/) const;
138 G4double CalculateChargeCollection(G4double, const G4ThreeVector&, const G4ThreeVector&, G4double /*Barret_PhiStart*/) const;
139 G4double CalculateChargeCollection1(G4double, const G4ThreeVector&, const G4ThreeVector&, G4double /*Barret_PhiStart*/) const;
140
141 // get power of gap in signal calculation
142 DoubleProperty m_GApower{this, "EMECGapPower", 1.4}; // used as const after init
143 inline G4double GApower() const { return m_GApower; };
144
145 // **************************************************************************
146 //Declaration of variables,functions for charge collection
147 //J.T
148 // **************************************************************************
149
150
151 //variables specific for wheel geometry
169
170
171 G4double m_ElectrodeFanHalfThickness = 0.0; // used as const after init
172 G4double m_FanEleThicknessOld = 0.0; // used as const after init
173 G4double m_FanEleFoldRadiusOld = 0.0; // used as const after init
174 G4double m_FanAbsThickness = 0.0; // used as const after init
175 G4double m_FanEleThickness = 0.0; // used as const after init
176 G4double m_WaveLength = 0.0; // used as const after init
177
178 G4double m_zsep12[44]{}; // used as const after initialization
179 G4double m_ziw[7]{}; // used as const after initialization
180 G4double m_zsep23[22]{}; // used as const after initialization
181
182
183 inline G4double ElectrodeFanHalfThickness() const { return m_ElectrodeFanHalfThickness; };
184 inline G4double FanEleThicknessOld() const { return m_FanEleThicknessOld; };
185 inline G4double FanEleFoldRadiusOld() const { return m_FanEleFoldRadiusOld; };
186 inline G4double FanAbsThickness() const { return m_FanAbsThickness; };
187 inline G4double FanEleThickness() const { return m_FanEleThickness; };
188 inline G4double WaveLength() const { return m_WaveLength; };
189
190 //variables specific for Efield calculation
191
192 G4String m_FieldMapVersion; // used as const after init
193
194 static const G4double s_GridSize;
195 static const G4double s_AverageGap;
196 static const G4double s_inv_AverageGap;
197
200 G4double* FieldMap; // [NumberOfRadialLayers][ZYWeight][MaxNofPoints];
201 G4double* MinZofLayer; //these are limits of the
202 G4double* MaxZofLayer; //area where the FieldMap can
203 G4double* MinYofLayer; //be used for interpolation
204 G4double* MaxYofLayer;
205 G4int* NofColofLayer; // a column is parallel to y
206 G4int* NofRowofLayer; // a row is parallel to z
208 G4int* pLayer;
209 };
210
213 G4double* RadiusOfLayers;
219 G4double GridShift;
220 };
221
222 Wheel_Efield_Map m_ChCollInner{},m_ChCollOuter{}; // used as const after init
223 Wheel_Efield_Map* m_ChCollWheelType = nullptr; // used as const after init
224
225 inline const Wheel_Efield_Map* ChCollWheelType() const { return m_ChCollWheelType; };
226
235
237 inline G4int Index(const Fold_Efield_Map* foldmap, G4int i, G4int j, G4int k ) const {
238 return foldmap->pLayer[i]+j*foldmap->NofPointsinLayer[i]+k;
239 };
240 void SetFoldArea(G4double, FoldArea & ) const;
241
242 StringProperty m_HVMapVersion{this, "EMECHVMap", "v02"}; // used only for initialization
243 BooleanProperty m_DB_HV{this, "EMECHVEnable", false};
244
245 static const G4double s_AverageHV;
246 static const G4double s_AverageEfield;
247 static const G4double s_AverageCurrent;
248
249 static const G4double s_LArTemperature_ECC0;//={88.15}; //K
250 static const G4double s_LArTemperature_ECC1;//={88.37};
251 static const G4double s_LArTemperature_ECC5;//={87.97};
252 static const G4double s_LArTemperature_av ;// ={88.16};
253
254 //Efield in [kv/cm], driftvelo in [mm/microsec], Temperature in [K]
255
256 inline static G4double IonReco(const G4double Efield) {
257 if(Efield<=0.000001){return 0.;}
258 if(Efield>2.) {return (1./(1. +0.36/Efield));}
259 return (1./(1.04+0.28/Efield));
260 }
261
262 inline static G4double DriftVelo(const G4double T, const G4double Efield) {
263 if( Efield <= 0.000001) {return 0.;}
264 return ( (-0.01481*(T-90.371)+1.)*
265 ( 0.141*Efield*log(1.+12.4/Efield)+
266 1.627*pow(Efield,0.317) )
267 -0.0075*(T-90.371)
268 );
269 }
270
271 // functions specific for geometry
272
273 void SetHalfWave(G4double, WheelGeometry &) const;
274 void GetPhiGap(const G4double *, WheelGeometry &) const;
275 void SetYlimitsofPhigapinWheel(G4double, G4double, const WheelGeometry & wg, G4double * Ylimits) const;
276 G4double YofSurface(G4double,G4double,G4double,const WheelGeometry &) const;
277 inline G4double YofNeutralFibre(G4double alpha,G4double rho, const WheelGeometry & wg) const {
278 return YofSurface(alpha,rho,0., wg);
279 }
280 G4double FoldingAngle(G4double) const;
281 G4double HalfLArGapSize(G4double, G4double) const;
282
283 // functions specific for charge coll.
284
285 void IniGeomforFieldMaps(void); // called only at init phase
286 void LoadFieldMaps(const G4String&); // called only at init phase
287 void PrepareFieldMap(Wheel_Efield_Map* ChCollWheelType); // called only at init phase
288 G4double GetCurrent(const G4double *,const G4double *,G4double, G4double Barret_PhiStart) const;
289 void TransformWheeltoFieldMap(const G4double *,G4double *, const WheelGeometry & wg, const FoldArea & fa) const;
290 void SetYlimitsofPhigapinFieldMap(G4int, const WheelGeometry & wg, G4double * Ylimits) const;
291 void TransFromBarrtoWheel(const G4double*, G4double PhiStartOfPhiDiv, G4double*) const;
292 G4double GetWeightfromFieldMap(G4int,G4double,G4double, const FoldArea & fa) const;
293 G4double HalfLArGapSizeOld(G4double) const;
294
295 // pick up surface_suppression_range
296 DoubleProperty m_CHC_Esr{this, "EMECRsr", 0.2*CLHEP::mm}; // used as const after init
297 inline G4double CHC_Esr() const { return m_CHC_Esr; };
298
299#ifdef DEBUG_CHCL // non thread-safe debug of charge collection
300 static const G4int s_CHCMaxPrint=00; // exists only if debug activated
301 static G4int s_CHCIprint; // exists only if debug activated
302 static G4double s_CHCEbad; // exists only if debug activated
303 static G4double s_CHCEtotal; // exists only if debug activated
304 static G4double s_CHCStotal; // exists only if debug activated
305#endif
306
307 private:
308 /* to be developed...
309 std::pair<double, double>DxToFans(Hep3Vector &p);
310 double XDistanceToTheNeutralFibre(const Hep3Vector& P) const;
311 */
312 G4double GetGapSize(const G4ThreeVector &p) const;
313
314 // public:
315 G4double distance_to_the_nearest_electrode(const G4ThreeVector &p) const;
316
317 UnsignedIntegerProperty m_solidtypeProp{this, "WheelType", 0, &EnergyCalculator::SolidTypeHandler};
319 IntegerProperty m_zside{this, "zSide", 1};
322 const LArWheelCalculator * lwc() const { return m_lwc; }
323
324 ServiceHandle<IGeoModelSvc> m_geoModel{this, "GeoModelSvc", "GeoModelSvc"};
325 ServiceHandle<IGeoDbTagSvc> m_geoDbTagSvc{this, "GeoDbTagSvc", "GeoDbTagSvc"};
326 StringProperty m_suffix{this, "EMECChMap", "v03"};
327
328 // Aug 2007 AMS, lost Aug 2008, restored May 2009
331
332 G4double GetCurrent1(const G4ThreeVector &, const G4ThreeVector &, G4double) const;
333
336
337 G4int _getIRlayer(G4double rforalpha) const;
338 G4int _getIRlayerA(G4double rforalpha) const;
339
340 G4double _interpolateCurrentSubStep(G4double rforalpha, G4int gapup, const G4double vmap[],
341 G4double tol, const FoldArea & fa, G4int & gaperr ) const;
342
343
344 G4double _interpolateCurrentSubStep1(G4double rforalpha, const G4double vmap[],
345 const G4ThreeVector & Pe, int side_dte, int Pe_fan,
346 const G4ThreeVector & Pa, int side_dta, int Pa_fan,
347 const FoldArea & fa, G4int & gaperr ) const;
348
349 G4double _AdjustedPhiOfPoint_Barrett(const G4ThreeVector& p, G4double PhiStartOfPhiDiv) const;
350
351 static inline G4double _normalizeAngle2Pi(G4double a) {
352 return ( a<0.) ?
353 a + CLHEP::twopi
354 :
355 (a >= CLHEP::twopi ? a - CLHEP::twopi: a);
356 }
357
358 G4double getPhiStartOfPhiDiv(const G4Step* step) const;
359
360 private:
361 std::unique_ptr<const HVHelper> m_HVHelper;
362 const static G4double s_GA_SubstepSize;
363 G4double DistanceToEtaLine(const G4ThreeVector &p, G4double eta) const;
364
365 struct geometry_t {
366 G4int zSide; // +- 3 for inner wheel, +- 2 for outer wheel, z determines sign
367 G4int sampling;
368 G4int region;
369 G4double etaScale; // 1/deta
370 G4double etaOffset; // set so that the range of etaBin starts at zero for each compartment
371 G4int maxEta; // the maximum value of etaBin in this compartment
372 G4int gapsPerBin; // number of phi gaps (in LArWheelSolid) for each cell bin.
373 G4int maxPhi; // the maximum value of phiBin in this compartment
374 };
375 static const geometry_t s_geometry[];
376
378 const G4ThreeVector& p, G4double PhiStartOfPhiDiv,
379 G4double &phi, G4int &compartment, G4int &eta_bin
380 ) const;
381
382 };
383 } // namespace EC
384} // namespace LArG4
385
386#endif // LArG4_EC_EnergyCalculator_H
float hitTime(const AFP_SIDSimHit &hit)
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
static Double_t a
constexpr int pow(int base, int exp) noexcept
LArCalculatorSvcImp(const std::string &name, ISvcLocator *pSvcLocator)
Gaudi::Property< double > m_OOTcut
G4bool Process_Barrett(const G4Step *, std::vector< LArHitData > &) const
G4double CalculateChargeCollection1(G4double, const G4ThreeVector &, const G4ThreeVector &, G4double) const
G4bool Process_Default(const G4Step *, std::vector< LArHitData > &) const
G4double distance_to_the_nearest_electrode(const G4ThreeVector &p) const
virtual StatusCode initialize() override final
G4double(EnergyCalculator::* m_GetGapSize_type)(const G4ThreeVector &p) const
static const G4double s_S3_Etalim[21]
Wheel_Efield_Map * m_ChCollWheelType
EnergyCalculator(const EnergyCalculator &)
G4double GetGapSize_Barrett(const G4ThreeVector &p) const
UnsignedIntegerProperty m_solidtypeProp
G4double GapAdjustment__sE(G4double, const G4ThreeVector &, const G4ThreeVector &, G4double) const
G4double DistanceToEtaLine(const G4ThreeVector &p, G4double eta) const
G4double YofSurface(G4double, G4double, G4double, const WheelGeometry &) const
static const G4double s_AverageEfield
G4int _getIRlayer(G4double rforalpha) const
G4double GapAdjustment_s(G4double, const G4ThreeVector &, const G4ThreeVector &, G4double) const
G4double GetWeightfromFieldMap(G4int, G4double, G4double, const FoldArea &fa) const
G4double dummy_correction_method(G4double e, const G4ThreeVector &, const G4ThreeVector &, G4double) const
void LoadFieldMaps(const G4String &)
void SetFoldArea(G4double, FoldArea &) const
static const G4double s_LArTemperature_ECC1
static const G4double s_AverageHV
void CreateArrays(Wheel_Efield_Map &, G4int)
static const G4double s_LArTemperature_ECC5
G4double _interpolateCurrentSubStep1(G4double rforalpha, const G4double vmap[], const G4ThreeVector &Pe, int side_dte, int Pe_fan, const G4ThreeVector &Pa, int side_dta, int Pa_fan, const FoldArea &fa, G4int &gaperr) const
static const G4double s_Zmeas_outer[2]
static const G4double s_ZmaxOfSignal
static const geometry_t s_geometry[]
G4double GetCurrent1(const G4ThreeVector &, const G4ThreeVector &, G4double) const
static const G4double s_GridSize
G4bool GetCompartment_Barrett(const G4ThreeVector &, G4double, G4double, G4double, G4int &, G4int &) const
static const G4double s_LArTemperature_ECC0
G4double HalfLArGapSizeOld(G4double) const
LArWheelCalculator * m_electrode_calculator
LArG4::LArWheelCalculator_t m_solidtype
G4bool GetVolumeIndex(const G4Step *, G4int &, G4int &) const
G4bool FindIdentifier_Default(const G4Step *, std::vector< LArHitData > &, G4ThreeVector &, G4ThreeVector &) const
void GetPhiGap(const G4double *, WheelGeometry &) const
G4double distance_to_the_nearest_electrode_Default(const G4ThreeVector &p, G4double) const
void CorrectionTypeHandler(Gaudi::Details::PropertyBase &)
ServiceHandle< ILArCalibCalculatorSvc > m_supportCalculator
static const G4double s_ColdCorrection
virtual G4bool isInTime(G4double hitTime) const override final
G4bool(EnergyCalculator::* m_Process_type)(const G4Step *, std::vector< LArHitData > &) const
G4double GapAdjustment_old(G4double, const G4ThreeVector &, const G4ThreeVector &, G4double) const
G4double GapAdjustment(G4double, const G4ThreeVector &, const G4ThreeVector &, G4double) const
G4bool FindIdentifier_Barrett(const G4Step *, G4double, std::vector< LArHitData > &, G4ThreeVector &, G4ThreeVector &) const
static const G4double s_AverageCurrent
G4double(EnergyCalculator::* m_ecorr_method)(G4double, const G4ThreeVector &, const G4ThreeVector &, G4double) const
G4double _AdjustedPhiOfPoint_Barrett(const G4ThreeVector &p, G4double PhiStartOfPhiDiv) const
const Wheel_Efield_Map * ChCollWheelType() const
virtual StatusCode finalize() override final
virtual G4bool Process(const G4Step *, std::vector< LArHitData > &) const override final
G4double(EnergyCalculator::* m_distance_to_the_nearest_electrode_type)(const G4ThreeVector &p, G4double) const
static const G4double s_KapGap
virtual G4float OOTcut() const override final
G4bool FindDMIdentifier_Barrett(const G4Step *step, std::vector< LArHitData > &) const
G4double FanEleFoldRadiusOld() const
static const G4double s_DistOfStartofCuFromBack
G4double GetGapSize(const G4ThreeVector &p) const
G4double GetGapSize_Default(const G4ThreeVector &p) const
static const G4double s_LongBarThickness
void SetYlimitsofPhigapinFieldMap(G4int, const WheelGeometry &wg, G4double *Ylimits) const
static const G4double s_StripWidth
G4double ElectrodeFanHalfThickness() const
G4bool GetBarrettePCE(const G4ThreeVector &p, G4double PhiStartOfPhiDiv, G4double &phi, G4int &compartment, G4int &eta_bin) const
G4double GetHV_Value(const G4ThreeVector &p, G4double PhiStartOfPhiDiv) const
static const G4double s_DistOfEndofCuFromBack
void SetHalfWave(G4double, WheelGeometry &) const
static const G4double s_AverageGap
G4double GetCurrent(const G4double *, const G4double *, G4double, G4double Barret_PhiStart) const
G4double CalculateChargeCollection(G4double, const G4ThreeVector &, const G4ThreeVector &, G4double) const
static const G4double s_Rmeas_outer[50]
G4double FoldingAngle(G4double) const
UnsignedIntegerProperty m_corrProp
G4double YofNeutralFibre(G4double alpha, G4double rho, const WheelGeometry &wg) const
static const G4double s_LArTemperature_av
EnergyCalculator(const std::string &name, ISvcLocator *pSvcLocator)
G4int Index(const Fold_Efield_Map *foldmap, G4int i, G4int j, G4int k) const
static const G4double s_EdgeWidth
void TransformWheeltoFieldMap(const G4double *, G4double *, const WheelGeometry &wg, const FoldArea &fa) const
G4double GapAdjustment_E(G4double, const G4ThreeVector &, const G4ThreeVector &, G4double) const
ServiceHandle< IGeoDbTagSvc > m_geoDbTagSvc
EnergyCalculator & operator=(const EnergyCalculator &)
static const G4double s_inv_AverageGap
void SolidTypeHandler(Gaudi::Details::PropertyBase &)
G4double HalfLArGapSize(G4double, G4double) const
void SetYlimitsofPhigapinWheel(G4double, G4double, const WheelGeometry &wg, G4double *Ylimits) const
void PrepareFieldMap(Wheel_Efield_Map *ChCollWheelType)
static G4double DriftVelo(const G4double T, const G4double Efield)
G4double FanEleThicknessOld() const
G4int GetPhiGap_Barrett(const G4ThreeVector &p, G4double PhiStartOfPhiDiv) const
static G4double IonReco(const G4double Efield)
G4int _getIRlayerA(G4double rforalpha) const
G4double distance_to_the_nearest_electrode_Barrett(const G4ThreeVector &p, G4double Barret_PhiStart) const
ServiceHandle< IGeoModelSvc > m_geoModel
G4double _interpolateCurrentSubStep(G4double rforalpha, G4int gapup, const G4double vmap[], G4double tol, const FoldArea &fa, G4int &gaperr) const
std::unique_ptr< const HVHelper > m_HVHelper
static const G4double s_GA_SubstepSize
void TransFromBarrtoWheel(const G4double *, G4double PhiStartOfPhiDiv, G4double *) const
const LArWheelCalculator * lwc() const
EnergyCorrection_t m_correction_type
G4double getPhiStartOfPhiDiv(const G4Step *step) const
static G4double _normalizeAngle2Pi(G4double a)
const LArWheelCalculator * elc() const
This class separates some of the geometry details of the LAr endcap.
EnergyCorrection_t
@ EMEC_ECOR_CHCL1