ATLAS Offline Software
Loading...
Searching...
No Matches
CryostatCalibrationLArCalculator.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// LArG4::BarrelCryostat::CalibrationLArCalculator
6// Prepared 06-Jul-2004 Bill Seligman
7// 16-Jul-2004 M.Leltchouk - extended for LAr::Barrel::MotherVolume
8
9// This class calculates the values needed for calibration hits in the
10// simulation.
11
12#undef DEBUG_HITS
13#define DEBUG_VOLUMES
14
16
18//#include "LArG4Code/CalibrationDefaultCalculator.h"
19
20#include "G4Step.hh"
21#include "G4StepPoint.hh"
22#include "G4VPhysicalVolume.hh"
23#include "G4ThreeVector.hh"
24#include "globals.hh"
25
26#include <atomic>
27#include <map>
28#include <algorithm>
29#include <cmath>
30#include <climits>
31
32namespace LArG4 {
33
34 namespace BarrelCryostat {
35
37 // Methods
39
40 CalibrationLArCalculator::CalibrationLArCalculator(const std::string& name, ISvcLocator *pSvcLocator)
41 : LArCalibCalculatorSvcImp(name, pSvcLocator)
42 {
43 }
44
46 // Get the default calculator (hopefully a temporary piece of code):
48
49 // Access source of detector parameters.
51 return StatusCode::SUCCESS;
52 }
53
54 G4bool CalibrationLArCalculator::Process(const G4Step* step, LArG4Identifier & identifier,
55 std::vector<G4double> & energies,
57 {
58 // Use the calculators to determine the energies and the
59 // identifier associated with this G4Step. Note that the
60 // default is to process both the energy and the ID.
61
63 {
64 m_energyCalculator.Energies( step, energies );
65 }
66 else
67 for (unsigned int i=0; i != 4; i++) energies.push_back( 0. );
68
69
70 identifier.clear();
71 if ( process == kEnergyAndID || process == kOnlyID )
72 {
73 // Calculate the identifier.
74
75 // Note:
76 // LArG4::BarrelCryostat::CryostatCalibrationCalculator uses
77 // a table-based approach to determine the identifier. The
78 // following code uses an "if-statement" approach.
79
80 // The fixed parameters (only a couple of which are readily
81 // accessible from the database):
82
83 constexpr double oneOverDeta = 10.; // 1/Deta = 1./0.1 = 10.
84 constexpr double oneOverDphi = 32./M_PI; // 1/Dphi
85
86 constexpr double rhoMinPresamplerMother = 1385.*CLHEP::mm;
87 constexpr double rhoMiddlePresampler = (1385.*CLHEP::mm + 1447.*CLHEP::mm)/2.;
88 // from PresParameterDef.icc
89 // rMinPresamplerMother = 1385*CLHEP::mm;
90 // rMaxPresamplerMother = 1447*CLHEP::mm - 0.001*CLHEP::mm;
91
92 constexpr double rhoAlignmentSafety = 10.*CLHEP::mm;
93 constexpr double rhoInFrontOfColdWall = rhoMinPresamplerMother - rhoAlignmentSafety;
94 static const double RIN_AC = m_parameters->GetValue("LArEMBRadiusInnerAccordion"); // 1500.024*CLHEP::mm; from ACCG
95 constexpr double RCUT12 = 1593.9*CLHEP::mm;
96 constexpr double RCUT23 = 1866.1*CLHEP::mm;
97 static const double ROUT_AC = m_parameters->GetValue("LArEMBRadiusOuterAccordion"); // 1960.*CLHEP::mm;
98 static const double rhoOuterAccordionWithSafety = ROUT_AC - rhoAlignmentSafety;
99 static const double LArEMBZmax = m_parameters->GetValue("LArEMBZmax"); // 3165.*CLHEP::mm
100
101 static const double zMaxAccordionWithSafety = LArEMBZmax - 10.*CLHEP::mm;
102
103
104 // Calculate the mid-point of the step, and the simple geometry variables.
105
106 G4StepPoint* pre_step_point = step->GetPreStepPoint();
107 G4StepPoint* post_step_point = step->GetPostStepPoint();
108
109 G4ThreeVector startPoint = pre_step_point->GetPosition();
110 G4ThreeVector endPoint = post_step_point->GetPosition();
111 G4ThreeVector p = (startPoint + endPoint) * 0.5;
112
113 G4double rho = p.perp();
114 G4double eta = fabs( p.pseudoRapidity() );
115 G4double phi = p.phi();
116 if ( phi < 0. ) phi += 2.*M_PI; // Normalize for phiBin calculation
117
118 // subdet = +/-4 "+" or " -" according to sign of Z in World coorinate
119 G4int subdet = ( p.z() > 0.) ? 4 : -4;
120 G4int phiBin = (int) ( phi * oneOverDphi );
121 if (phiBin>63) phiBin=0;
122
123 // Initialize identifier variables with (invalid) default
124 // values (INT_MIN is defined in <climits>).
125 G4int type = INT_MIN;
126 G4int sampling = INT_MIN;
127 G4int region = INT_MIN;
128 G4int etaBin = INT_MIN;
129
130
131 if(std::fabs(p.z())< 10 && eta < 0.1 && rho > rhoMinPresamplerMother && rho < ROUT_AC )
132 {
133 // type = 2 dead materials in "internal" cracks:
134 // real crack for nominal geometry at this level is +-3mm
135 // add some margin to allow for misalignement
136
137 type = 2;
138
139 // sampling = 0,1,2,3 (follow the corresponding calorimeter structure)
140 if ( rho < RIN_AC ) { sampling = 0; }
141 else if ( rho < RCUT12 ) { sampling = 1; }
142 else if ( rho < RCUT23 ) { sampling = 2; }
143 else { sampling = 3; }
144
145 // region = 0 - between two halves of EMB,
146 region = 0;
147 etaBin = 0;
148 }
149 else
150 {
151 // type = 1 dead materials outside accordion and active presampler layers
152 type = 1;
153
154 if( eta < 1.5 )
155 {
156 if ( rho < rhoInFrontOfColdWall ) // and in LAr::Barrel::MotherVolume
157
158 // The region is defined correctly for a radial mis-alignments
159 // up to 10mm because ColdWallThickness >= 20mm and
160 // rhoAlignmentSafety = 10.*mm;
161 // The E-deposit in the Cold Wall itself is handled by
162 // another calculator.
163 {
164 sampling = 1;
165 region = 0;
166 etaBin = (int) ( eta * oneOverDeta );
167 }
168 else if ( rho < rhoMiddlePresampler )
169 // The region is defined correctly for a radial mis-alignments
170 // up to (1447.*mm - 1385.*mm)/2
171 // The E-deposit in the Presampler itself is handled by
172 // other calculator.
173 {
174 sampling = 1;
175 region = 2;
176 etaBin = (int) ( eta * oneOverDeta );
177 }
178 else if ( rho < rhoOuterAccordionWithSafety &&
179 fabs( p.z()) < zMaxAccordionWithSafety )
180 // The region is defined correctly for a radial mis-alignments
181 // up to 10mm because rhoAlignmentSafety = 10.*mm;
182 // The E-deposit in the Barrel Accordion itself is handled by
183 // other calculator.
184 {
185 sampling = 1;
186 region = 3;
187 etaBin = (int) ( eta * oneOverDeta );
188 }
189 else // rho > ROUT_AC or |z| >= zMaxAccordionWithSafety
190 {
191 sampling = 2;
192
193 if( eta < 1.0 )
194 {
195 region = 0;
196 etaBin = (int) ( eta * oneOverDeta );
197 }
198 else // 1.0 <= eta < 1.5
199 {
200 region = 2;
201 etaBin = (int) ( (eta-1.) * oneOverDeta );
202 }
203 }
204 }
205 else if ( eta < 1.6 )
206 {
207 sampling = 1;
208 region = 4;
209 etaBin = (int) ( (eta-1.5) * oneOverDeta );
210 }
211 else if ( eta < 1.8 )
212 {
213 sampling = 1;
214 region = 5;
215 etaBin = (int) ( (eta-1.6) * oneOverDeta );
216 }
217 else if ( eta < 3.2 )
218 {
219 sampling = 1;
220 region = 6;
221 etaBin = (int) ( (eta-1.8) * oneOverDeta );
222 }
223 }
224
225 // This is a "quick fix" for a complex issue: We're still
226 // developing code for the cryostat. What if, somehow, we have
227 // a G4Step in a LAr volume that isn't handled by the above
228 // code? Answer: Use the default calibration calculator (the
229 // same one used for volumes without sensitive detectors) to get
230 // the identifier.
231
232 if ( type == INT_MIN ||
233 region == INT_MIN ||
234 sampling == INT_MIN ||
235 etaBin == INT_MIN ||
236 phiBin < 0 )
237 {
238#if defined (DEBUG_VOLUMES) || defined (DEBUG_HITS)
239 constexpr G4int messageMax = 10;
240 static std::atomic<G4int> messageCount = 0;
241 if ( messageCount++ < messageMax )
242 {
243 std::cout << "LArG4::BarrelCryostat::CalibrationLArCalculator::Process"
244 << " (error " << messageCount << " of " << messageMax << " max displayed)"
245 << std::endl
246 << " G4Step in LAr at unexpected place: (x,y,z) [mm] = ("
247 << p.x()/CLHEP::mm << ","
248 << p.y()/CLHEP::mm << ","
249 << p.z()/CLHEP::mm
250 << "); eta=" << eta
251 << ", phi=" << phi << std::endl
252 << " using default calculator"
253 << std::endl;
254 }
255#endif
256 m_defaultCalculator->Process(step, identifier, energies, process);
257 }
258 else
259 {
260 // Append the cell ID to the (empty) identifier.
261 identifier << 10 // Calorimeter
262 << subdet // LAr +/-4 where "+" or " -" according to
263 // the sign of Z in World coorinate
264 << type
265 << sampling
266 << region
267 << etaBin
268 << phiBin;
269 }
270 }
271
272#ifdef DEBUG_HITS
273 G4double energy = accumulate(energies.begin(),energies.end(),0.);
274 std::cout << "LArG4::BarrelCryostat::CalibrationLArCalculator::Process"
275 << " ID=" << std::string(identifier)
276 << " energy=" << energy
277 << " energies=(" << energies[0]
278 << "," << energies[1]
279 << "," << energies[2]
280 << "," << energies[3] << ")"
281 << std::endl;
282#endif
283
284 // Check for bad result.
285 if ( identifier == LArG4Identifier() )
286 return false;
287
288 return true;
289 }
290
291 } // namespace BarrelCryostat
292
293} // namespace LArG4
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
bool accumulate(AccumulateMap &map, std::vector< module_t > const &modules, FPGATrackSimMatrixAccumulator const &acc)
Accumulates an accumulator (e.g.
LArCalibCalculatorSvcImp(const std::string &name, ISvcLocator *pSvcLocator)
virtual G4bool Process(const G4Step *step, LArG4Identifier &identifier, std::vector< G4double > &energies, const eCalculatorProcessing process=kEnergyAndID) const override final
CalibrationLArCalculator(const std::string &name, ISvcLocator *pSvcLocator)
static const VDetectorParameters * GetInstance()
const std::string process
eCalculatorProcessing
@ kOnlyEnergy
@ kEnergyAndID