ATLAS Offline Software
Loading...
Searching...
No Matches
LArG4::BarrelCryostat::CalibrationLArCalculator Class Reference

#include <CryostatCalibrationLArCalculator.h>

Inheritance diagram for LArG4::BarrelCryostat::CalibrationLArCalculator:
Collaboration diagram for LArG4::BarrelCryostat::CalibrationLArCalculator:

Public Member Functions

 CalibrationLArCalculator (const std::string &name, ISvcLocator *pSvcLocator)
StatusCode initialize () override final
virtual ~CalibrationLArCalculator ()=default
virtual G4bool Process (const G4Step *step, LArG4Identifier &identifier, LArG4Identifier &identifier_sr, std::vector< double > &energies, const LArG4::eCalculatorProcessing process) const override final

Private Attributes

CaloG4::SimulationEnergies m_energyCalculator {}
const LArVG4DetectorParametersm_parameters {nullptr}
ServiceHandle< ILArCalibCalculatorSvcm_defaultCalculator {this, "DefaultCalculator", "CalibrationDefaultCalculator"}

Detailed Description

Constructor & Destructor Documentation

◆ CalibrationLArCalculator()

LArG4::BarrelCryostat::CalibrationLArCalculator::CalibrationLArCalculator ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 40 of file CryostatCalibrationLArCalculator.cxx.

41 : LArCalibCalculatorSvcImp(name, pSvcLocator)
42 {
43 }
LArCalibCalculatorSvcImp(const std::string &name, ISvcLocator *pSvcLocator)

◆ ~CalibrationLArCalculator()

virtual LArG4::BarrelCryostat::CalibrationLArCalculator::~CalibrationLArCalculator ( )
virtualdefault

Member Function Documentation

◆ initialize()

StatusCode LArG4::BarrelCryostat::CalibrationLArCalculator::initialize ( )
finaloverride

Definition at line 45 of file CryostatCalibrationLArCalculator.cxx.

45 {
46 // Get the default calculator (hopefully a temporary piece of code):
48
49 // Access source of detector parameters.
51 return StatusCode::SUCCESS;
52 }
#define ATH_CHECK
Evaluate an expression and check for errors.
static const VDetectorParameters * GetInstance()

◆ Process()

G4bool LArG4::BarrelCryostat::CalibrationLArCalculator::Process ( const G4Step * step,
LArG4Identifier & identifier,
LArG4Identifier & identifier_sr,
std::vector< double > & energies,
const LArG4::eCalculatorProcessing process ) const
finaloverridevirtual

Definition at line 54 of file CryostatCalibrationLArCalculator.cxx.

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

Member Data Documentation

◆ m_defaultCalculator

ServiceHandle<ILArCalibCalculatorSvc> LArG4::BarrelCryostat::CalibrationLArCalculator::m_defaultCalculator {this, "DefaultCalculator", "CalibrationDefaultCalculator"}
private

Definition at line 77 of file LArG4Barrel/src/CryostatCalibrationLArCalculator.h.

77{this, "DefaultCalculator", "CalibrationDefaultCalculator"};

◆ m_energyCalculator

CaloG4::SimulationEnergies LArG4::BarrelCryostat::CalibrationLArCalculator::m_energyCalculator {}
private

Definition at line 72 of file LArG4Barrel/src/CryostatCalibrationLArCalculator.h.

72{};

◆ m_parameters

const LArVG4DetectorParameters* LArG4::BarrelCryostat::CalibrationLArCalculator::m_parameters {nullptr}
private

Definition at line 75 of file LArG4Barrel/src/CryostatCalibrationLArCalculator.h.

75{nullptr};

The documentation for this class was generated from the following files: