ATLAS Offline Software
Loading...
Searching...
No Matches
LArWheelSliceSolidInit.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <cassert>
6#include <stdexcept>
7#include <iostream>
8
9#ifndef PORTABLE_LAR_SHAPE
13
18
19
20#include "GaudiKernel/Bootstrap.h"
21#include "GaudiKernel/ISvcLocator.h"
22#include "GaudiKernel/IMessageSvc.h"
23#endif
24
26#include "CLHEP/Units/PhysicalConstants.h"
27#ifndef PORTABLE_LAR_SHAPE
29#endif
30#include "G4GeometryTolerance.hh"
31
33#include "LArWheelSliceSolid.h"
34#include "G4ShiftedCone.h"
35
36#ifdef DEBUG_LARWHEELSLICESOLID
37G4int LArWheelSliceSolid::Verbose = 0;
38#endif
39
40// these are internal technical constants, should not go in DB
41const unsigned int LArWheelSliceSolid::s_IterationsLimit = 50; // That's enough even for 10e-15 IterationPrecision
42const G4double LArWheelSliceSolid::s_Tolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance() / 2;
43const G4double LArWheelSliceSolid::s_AngularTolerance = G4GeometryTolerance::GetInstance()->GetAngularTolerance() / 2;
44const G4double LArWheelSliceSolid::s_IterationPrecision = 0.001*CLHEP::mm;
45const G4double LArWheelSliceSolid::s_IterationPrecision2 = s_IterationPrecision * s_IterationPrecision;
46
48 const G4String& name,
49 pos_t pos, type_t type, size_t slice,
50 G4int zside,
51 const LArWheelCalculator *calc,
52 const EMECData *emecData
53) : G4VSolid(name)
54#ifndef PORTABLE_LAR_SHAPE
55 , AthMessaging("LArWSS")
56#endif
57 , m_Pos(pos), m_Type(type), m_Calculator(calc)
58{
59 createSolid(name, zside, slice, emecData);
60}
61
62LArWheelSliceSolid::LArWheelSliceSolid(const G4String& name, const EMECData *emecData)
63 : G4VSolid(name)
64#ifndef PORTABLE_LAR_SHAPE
65 , AthMessaging("LArWSS")
66#endif
67 , m_Calculator(0)
68{
69 if(name.find("::Inner") != G4String::npos) m_Pos = Inner;
70 else if(name.find("::Outer") != G4String::npos) m_Pos = Outer;
71 else G4Exception(
72 "LArWheelSliceSolid", "NoPos", FatalException,
73 (std::string("Constructor: can't get Inner/Outer from ") + name).c_str()
74 );
75
76 if(name.find("::Absorber") != G4String::npos) m_Type = Absorber;
77 else if(name.find("::Electrode") != G4String::npos) m_Type = Electrode;
78 else if(name.find("::Glue") != G4String::npos) m_Type = Glue;
79 else if(name.find("::Lead") != G4String::npos) m_Type = Lead;
80 else G4Exception(
81 "LArWheelSliceSolid", "NoType", FatalException,
82 (std::string("Constructor: can't get Type from ") + name).c_str()
83 );
84
85 size_t s = name.find("Slice");
86 size_t slice = 0;
87 if(G4String::npos != s){
88 slice = (size_t) atoi(name.substr(s + 5, 2).c_str());
89 } else G4Exception(
90 "LArWheelSliceSolid", "NoSlice", FatalException,
91 (std::string("Constructor: can't get Slice from ") + name).c_str()
92 );
93
94 G4int zside = 0;
95 if(name.find("::Pos::") != G4String::npos) zside = 1;
96 else if(name.find("::Neg::") != G4String::npos) zside = -1;
97 else G4Exception(
98 "LArWheelSliceSolid", "NoSide", FatalException,
99 (std::string("Constructor: can't get zSide from ") + name).c_str()
100 );
101 createSolid(name, zside, slice, emecData);
102}
103
104void LArWheelSliceSolid::createSolid(const G4String& name, G4int zside, size_t slice, const EMECData *emecData)
105{
106#ifndef PORTABLE_LAR_SHAPE
108#endif
109
111 switch(m_Pos){
112 case Inner:
113 switch(m_Type){
114 case Absorber: calc_type = LArG4::InnerAbsorberWheel; break;
115 case Electrode: calc_type = LArG4::InnerElectrodWheel; break;
116 case Glue: calc_type = LArG4::InnerGlueWheel; break;
117 case Lead: calc_type = LArG4::InnerLeadWheel; break;
118 }
119 break;
120 case Outer:
121 switch(m_Type){
122 case Absorber: calc_type = LArG4::OuterAbsorberWheel; break;
123 case Electrode: calc_type = LArG4::OuterElectrodWheel; break;
124 case Glue: calc_type = LArG4::OuterGlueWheel; break;
125 case Lead: calc_type = LArG4::OuterLeadWheel; break;
126 }
127 break;
128 }
129 if(m_Calculator == 0) {
130 m_Calculator = new LArWheelCalculator(*emecData, calc_type, zside);
131
132 }
133 else if(m_Calculator->type() != calc_type){
134 G4Exception(
135 "LArWheelSliceSolid", "WrongCalculatorType", FatalException,
136 "Constructor: external LArWheelCalculator of wrong type provided"
137 );
138 }
139
140 const G4String bs_name = name + "-Bounding";
141
142#ifdef DEBUG_LARWHEELSLICESOLID
143 const char *venv = getenv("LARWHEELSLICESOLID_VERBOSE");
144 if(venv) Verbose = atoi(venv);
145 static bool first = true;
146 if(first){
147 std::cout << "The LArWheelSliceSolid build "
148 << __DATE__ << " " << __TIME__ << std::endl
149 << "LArWheelSliceSolid verbosity level is "
150 << Verbose << std::endl;
151 first = false;
152 }
153#endif
154
158
159 switch(m_Pos){
160 case Inner: inner_solid_init(bs_name, slice); break;
161 case Outer: outer_solid_init(bs_name, slice); break;
162 }
163#ifndef PORTABLE_LAR_SHAPE
164 ATH_MSG_DEBUG(m_BoundingShape->GetName() + " is the m_BoundingShape");
165
166 init_tests();
167 test(); // activated by env. variable
168 clean_tests();
169
170 ATH_MSG_DEBUG("slice " << m_Pos << " " << m_Type
171 << " " << slice << " initialized" << endmsg);
172
173#endif
174
175 #ifdef DEBUG_LARWHEELSLICESOLID
176 std::cout << "LArWSS(" << m_Pos << ", " << m_Type
177 << "): slice " << slice << ", Zmin = " << m_Zmin
178 << ", Zmax = " << m_Zmax << std::endl
179 << GetName() << std::endl;
180#endif
181}
182
183inline void LArWheelSliceSolid::check_slice(size_t size, size_t slice) const
184{
185 if(size <= slice + 1){
186 G4Exception(
187 "LArWheelSliceSolid", "WrongSlice", FatalException,
188 std::string("Constructor: Slice number too big: " + GetName()).c_str()
189 );
190 }
191}
192
193void LArWheelSliceSolid::inner_solid_init(const G4String &bs_name, size_t slice)
194{
195 std::array<G4double,2> zPlane{}, rInner{}, rOuter{};
196 zPlane[1] = GetCalculator()->GetWheelThickness();
199
200 std::vector<G4double> zsect;
201 fill_zsect(zsect);
202 check_slice(zsect.size(), slice);
203
204 m_Zmin = zsect[slice]; m_Zmax = zsect[slice + 1];
205 m_Rmin = rInner[0]; m_Rmax = rOuter[1];
206 const G4double ainn = (rInner[1] - rInner[0]) / (zPlane[1] - zPlane[0]);
207 const G4double aout = (rOuter[1] - rOuter[0]) / (zPlane[1] - zPlane[0]);
208 const G4double R1inn = ainn * (m_Zmin - zPlane[0]) + rInner[0];
209 const G4double R1out = aout * (m_Zmin - zPlane[0]) + rOuter[0];
210 const G4double R2inn = ainn * (m_Zmax - zPlane[0]) + rInner[0];
211 const G4double R2out = aout * (m_Zmax - zPlane[0]) + rOuter[0];
212
214 bs_name + "Cone", m_Zmin, m_Zmax, R1inn, R1out, R2inn, R2out
215 );
216
217 const G4double FanPhiAmplitude = 0.065; // internal technical constant, should not go in DB
218 const G4double phi_min = CLHEP::halfpi
219 - FanPhiAmplitude
221 m_Xmax = rOuter[1]*cos(phi_min);
222 m_Xmin = -m_Xmax;
223
224 m_Ymin = m_Rmin * 0.9;
225}
226
227void LArWheelSliceSolid::outer_solid_init(const G4String &bs_name, size_t slice)
228{
229 std::array<G4double,3> zPlane{}, rInner{}, rOuter{};
230 zPlane[1] = GetCalculator()->GetWheelInnerRadius(rInner);
231 zPlane[2] = GetCalculator()->GetWheelThickness();
233
234 std::vector<G4double> zsect;
235 fill_zsect(zsect, zPlane[1]);
236 check_slice(zsect.size(), slice);
237
238 m_Zmin = zsect[slice]; m_Zmax = zsect[slice + 1];
239 m_Rmin = rInner[0]; m_Rmax = rOuter[1];
240 const size_t i = m_Zmax > zPlane[1]? 1: 0;
241 const G4double dz = zPlane[i + 1] - zPlane[i];
242 const G4double ainn = (rInner[i + 1] - rInner[i]) / dz;
243 const G4double aout = (rOuter[i + 1] - rOuter[i]) / dz;
244 const G4double R1inn = ainn * (m_Zmin - zPlane[i]) + rInner[i];
245 const G4double R1out = aout * (m_Zmin - zPlane[i]) + rOuter[i];
246 const G4double R2inn = ainn * (m_Zmax - zPlane[i]) + rInner[i];
247 const G4double R2out = aout * (m_Zmax - zPlane[i]) + rOuter[i];
248
250 bs_name + "Cone", m_Zmin, m_Zmax, R1inn, R1out, R2inn, R2out
251 );
252
253 const G4double FanPhiAmplitude = 0.02; // internal technical constant, should not go in DB
254 const G4double phi_min = CLHEP::halfpi
255 - FanPhiAmplitude
257 m_Xmax = rOuter[1]*cos(phi_min);
258 m_Xmin = -m_Xmax;
259
260 m_Ymin = m_Rmin * 0.9;
261}
262
263void LArWheelSliceSolid::fill_zsect(std::vector<G4double> &zsect, G4double zMid) const
264{
265 const G4double half_wave_length = GetCalculator()->GetHalfWaveLength();
266 const G4double sss = GetCalculator()->GetStraightStartSection();
267 const G4int num_fs = GetCalculator()->GetNumberOfHalfWaves() + 1;
268 const G4double wheel_thickness = GetCalculator()->GetWheelThickness();
269
270 zsect.push_back(0.);
271 zsect.push_back(sss + half_wave_length * 0.25);
272 for(G4int i = 2; i < num_fs; i ++){
273 const G4double zi = half_wave_length * (i - 1) + sss;
274 if(m_Pos == Outer && zsect.back() < zMid && zi > zMid){
275 zsect.push_back(zMid);
276 }
277 zsect.push_back(zi);
278 }
279 zsect.push_back(wheel_thickness - sss - half_wave_length * 0.25);
280 zsect.push_back(wheel_thickness);
281}
#define endmsg
#define ATH_MSG_DEBUG(x)
Definition of the abstract IRDBAccessSvc interface.
Definition of the abstract IRDBRecord interface.
Definition of the abstract IRDBRecordset interface.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
This class separates some of the geometry details of the LAr endcap.
double GetFanStepOnPhi() const
int GetNumberOfHalfWaves() const
double GetFanHalfThickness(LArG4::LArWheelCalculator_t) const
double GetHalfWaveLength() const
double GetStraightStartSection() const
double GetWheelThickness() const
double GetWheelInnerRadius(std::array< double, 2 > &rInner) const
void GetWheelOuterRadius(std::array< double, 2 > &rOuter) const
static const G4double s_IterationPrecision
const LArWheelCalculator * m_Calculator
void check_slice(size_t size, size_t slice) const
void fill_zsect(std::vector< G4double > &, G4double zMid=0.) const
static const unsigned int s_IterationsLimit
LArWheelSliceSolid(const G4String &name, pos_t pos, type_t type, size_t slice, G4int zside=1, const LArWheelCalculator *calc=0, const EMECData *emecData=0)
void inner_solid_init(const G4String &, size_t slice)
const LArWheelCalculator * GetCalculator(void) const
static const G4double s_Tolerance
static const G4double s_IterationPrecision2
void createSolid(const G4String &name, G4int zside, size_t slice, const EMECData *emecData)
static const G4double s_AngularTolerance
void outer_solid_init(const G4String &, size_t slice)